xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision d9558ea91a0ca1fa4633d87e21172ffbe9a39015)
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*/
6b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
72e68589bSMark F. Adams 
82e68589bSMark F. Adams #include <petscblaslapack.h>
92e68589bSMark F. Adams 
102e68589bSMark F. Adams typedef struct {
11c8b0795cSMark F. Adams   PetscInt  nsmooths;
12c8b0795cSMark F. Adams   PetscBool sym_graph;
13ef4ad70eSMark F. Adams   PetscBool square_graph;
142e68589bSMark F. Adams } PC_GAMG_AGG;
152e68589bSMark F. Adams 
162e68589bSMark F. Adams #undef __FUNCT__
172e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths"
182e68589bSMark F. Adams /*@
192e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
202e68589bSMark F. Adams 
212e68589bSMark F. Adams    Not Collective on PC
222e68589bSMark F. Adams 
232e68589bSMark F. Adams    Input Parameters:
242e68589bSMark F. Adams .  pc - the preconditioner context
252e68589bSMark F. Adams 
262e68589bSMark F. Adams    Options Database Key:
272e68589bSMark F. Adams .  -pc_gamg_agg_nsmooths
282e68589bSMark F. Adams 
292e68589bSMark F. Adams    Level: intermediate
302e68589bSMark F. Adams 
312e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
322e68589bSMark F. Adams 
332e68589bSMark F. Adams .seealso: ()
342e68589bSMark F. Adams @*/
352e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
362e68589bSMark F. Adams {
372e68589bSMark F. Adams   PetscErrorCode ierr;
382e68589bSMark F. Adams 
392e68589bSMark F. Adams   PetscFunctionBegin;
402e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
412e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
422e68589bSMark F. Adams   PetscFunctionReturn(0);
432e68589bSMark F. Adams }
442e68589bSMark F. Adams 
452e68589bSMark F. Adams #undef __FUNCT__
462e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG"
472e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n)
482e68589bSMark F. Adams {
492e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
502e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
51c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
522e68589bSMark F. Adams 
532e68589bSMark F. Adams   PetscFunctionBegin;
54c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
55c8b0795cSMark F. Adams   PetscFunctionReturn(0);
56c8b0795cSMark F. Adams }
57c8b0795cSMark F. Adams 
58c8b0795cSMark F. Adams #undef __FUNCT__
59c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph"
60c8b0795cSMark F. Adams /*@
61c8b0795cSMark F. Adams    PCGAMGSetSymGraph -
62c8b0795cSMark F. Adams 
63c8b0795cSMark F. Adams    Not Collective on PC
64c8b0795cSMark F. Adams 
65c8b0795cSMark F. Adams    Input Parameters:
66c8b0795cSMark F. Adams .  pc - the preconditioner context
67c8b0795cSMark F. Adams 
68c8b0795cSMark F. Adams    Options Database Key:
69c8b0795cSMark F. Adams .  -pc_gamg_sym_graph
70c8b0795cSMark F. Adams 
71c8b0795cSMark F. Adams    Level: intermediate
72c8b0795cSMark F. Adams 
73c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
74c8b0795cSMark F. Adams 
75c8b0795cSMark F. Adams .seealso: ()
76c8b0795cSMark F. Adams @*/
77c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
78c8b0795cSMark F. Adams {
79c8b0795cSMark F. Adams   PetscErrorCode ierr;
80c8b0795cSMark F. Adams 
81c8b0795cSMark F. Adams   PetscFunctionBegin;
82c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
83c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
84c8b0795cSMark F. Adams   PetscFunctionReturn(0);
85c8b0795cSMark F. Adams }
86c8b0795cSMark F. Adams 
87c8b0795cSMark F. Adams #undef __FUNCT__
88c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph_GAMG"
89c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n)
90c8b0795cSMark F. Adams {
91c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
92c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
93c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
94c8b0795cSMark F. Adams 
95c8b0795cSMark F. Adams   PetscFunctionBegin;
96c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
972e68589bSMark F. Adams   PetscFunctionReturn(0);
982e68589bSMark F. Adams }
992e68589bSMark F. Adams 
100ef4ad70eSMark F. Adams #undef __FUNCT__
101ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph"
102ef4ad70eSMark F. Adams /*@
103ef4ad70eSMark F. Adams    PCGAMGSetSquareGraph -
104ef4ad70eSMark F. Adams 
105ef4ad70eSMark F. Adams    Not Collective on PC
106ef4ad70eSMark F. Adams 
107ef4ad70eSMark F. Adams    Input Parameters:
108ef4ad70eSMark F. Adams .  pc - the preconditioner context
109ef4ad70eSMark F. Adams 
110ef4ad70eSMark F. Adams    Options Database Key:
111ef4ad70eSMark F. Adams .  -pc_gamg_square_graph
112ef4ad70eSMark F. Adams 
113ef4ad70eSMark F. Adams    Level: intermediate
114ef4ad70eSMark F. Adams 
115ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
116ef4ad70eSMark F. Adams 
117ef4ad70eSMark F. Adams .seealso: ()
118ef4ad70eSMark F. Adams @*/
119ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n)
120ef4ad70eSMark F. Adams {
121ef4ad70eSMark F. Adams   PetscErrorCode ierr;
122ef4ad70eSMark F. Adams 
123ef4ad70eSMark F. Adams   PetscFunctionBegin;
124ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125ef4ad70eSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
126ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
127ef4ad70eSMark F. Adams }
128ef4ad70eSMark F. Adams 
129ef4ad70eSMark F. Adams #undef __FUNCT__
130ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG"
131ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n)
132ef4ad70eSMark F. Adams {
133ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
134ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
135ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
136ef4ad70eSMark F. Adams 
137ef4ad70eSMark F. Adams   PetscFunctionBegin;
138ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
139ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
140ef4ad70eSMark F. Adams }
141ef4ad70eSMark F. Adams 
1422e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1432e68589bSMark F. Adams /*
1442e68589bSMark F. Adams    PCSetFromOptions_GAMG_AGG
1452e68589bSMark F. Adams 
1462e68589bSMark F. Adams   Input Parameter:
1472e68589bSMark F. Adams    . pc -
1482e68589bSMark F. Adams */
1492e68589bSMark F. Adams #undef __FUNCT__
1502e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
1518c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptions *PetscOptionsObject,PC pc)
1522e68589bSMark F. Adams {
1532e68589bSMark F. Adams   PetscErrorCode ierr;
1542e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1552e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
156c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1572e68589bSMark F. Adams 
1582e68589bSMark F. Adams   PetscFunctionBegin;
159e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1602e68589bSMark F. Adams   {
1612e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
162d3042614SMark Adams     pc_gamg_agg->nsmooths = 1;
1632fa5cd67SKarl Rupp 
1648afaa268SBarry 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);
165c8b0795cSMark F. Adams 
166c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
167c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
1682fa5cd67SKarl Rupp 
1698afaa268SBarry 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);
170ef4ad70eSMark F. Adams 
171ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
172ef4ad70eSMark F. Adams     pc_gamg_agg->square_graph = PETSC_TRUE;
1732fa5cd67SKarl Rupp 
1748afaa268SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_square_graph","For faster coarsening and lower coarse grid complexity","PCGAMGSetSquareGraph",pc_gamg_agg->square_graph,&pc_gamg_agg->square_graph,NULL);CHKERRQ(ierr);
1752e68589bSMark F. Adams   }
1762e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1772e68589bSMark F. Adams   PetscFunctionReturn(0);
1782e68589bSMark F. Adams }
1792e68589bSMark F. Adams 
1802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1812e68589bSMark F. Adams /*
1822e68589bSMark F. Adams    PCDestroy_AGG
1832e68589bSMark F. Adams 
1842e68589bSMark F. Adams   Input Parameter:
1852e68589bSMark F. Adams    . pc -
1862e68589bSMark F. Adams */
1872e68589bSMark F. Adams #undef __FUNCT__
1889b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG"
1899b8ffb57SJed Brown PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1902e68589bSMark F. Adams {
1912e68589bSMark F. Adams   PetscErrorCode ierr;
1922e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1932e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1942e68589bSMark F. Adams 
1952e68589bSMark F. Adams   PetscFunctionBegin;
1969b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1973ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
1982e68589bSMark F. Adams   PetscFunctionReturn(0);
1992e68589bSMark F. Adams }
2002e68589bSMark F. Adams 
2012e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2022e68589bSMark F. Adams /*
2032e68589bSMark F. Adams    PCSetCoordinates_AGG
204302f38e8SMark F. Adams      - collective
2052e68589bSMark F. Adams 
2062e68589bSMark F. Adams    Input Parameter:
2072e68589bSMark F. Adams    . pc - the preconditioner context
208a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
209302f38e8SMark F. Adams    . a_nloc - number of vertices local
210302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2112e68589bSMark F. Adams */
2121e6b0712SBarry Smith 
2132e68589bSMark F. Adams #undef __FUNCT__
2142e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2151e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2162e68589bSMark F. Adams {
2172e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2182e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2192e68589bSMark F. Adams   PetscErrorCode ierr;
22069344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
221a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2222e68589bSMark F. Adams 
2232e68589bSMark F. Adams   PetscFunctionBegin;
224a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
225a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
226302f38e8SMark F. Adams   nloc = a_nloc;
2272e68589bSMark F. Adams 
2282e68589bSMark F. Adams   /* SA: null space vectors */
22969344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
23069344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
231a2f3521dSMark F. Adams   else if (coords) {
232c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
23369344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
23469344418SMark F. Adams     if (ndm != ndf) {
235c666adbfSMark 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);
236a2f3521dSMark F. Adams     }
237806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
23871959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
23971959b99SBarry 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);
240c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2412e68589bSMark F. Adams 
2422e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2432e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2442e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
245854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2462e68589bSMark F. Adams   }
2472e68589bSMark F. Adams   /* copy data in - column oriented */
2482e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
24969344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
25069344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
251c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2522e68589bSMark F. Adams     else {
25369344418SMark F. Adams       /* translational modes */
2542fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2552fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
25669344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2572e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2582fa5cd67SKarl Rupp         }
2592fa5cd67SKarl Rupp       }
26069344418SMark F. Adams 
26169344418SMark F. Adams       /* rotational modes */
2622e68589bSMark F. Adams       if (coords) {
26369344418SMark F. Adams         if (ndm == 2) {
2642e68589bSMark F. Adams           data   += 2*M;
2652e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2662e68589bSMark F. Adams           data[1] =  coords[2*kk];
267806fa848SBarry Smith         } else {
2682e68589bSMark F. Adams           data   += 3*M;
2692e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2702e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2712e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2722e68589bSMark F. Adams         }
2732e68589bSMark F. Adams       }
2742e68589bSMark F. Adams     }
2752e68589bSMark F. Adams   }
2762e68589bSMark F. Adams 
2772e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2782e68589bSMark F. Adams   PetscFunctionReturn(0);
2792e68589bSMark F. Adams }
2802e68589bSMark F. Adams 
281b43b03e9SMark F. Adams typedef PetscInt NState;
282b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
283b43b03e9SMark F. Adams static const NState DELETED =-1;
284b43b03e9SMark F. Adams static const NState REMOVED =-3;
285b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
286b43b03e9SMark F. Adams 
287c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
288c8b0795cSMark F. Adams /*
289b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
290b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
291c8b0795cSMark F. Adams 
292c8b0795cSMark F. Adams    Input Parameter:
2932adfe9d3SBarry Smith    . Gmat_2 - glabal matrix of graph (data not defined)   base (squared) graph
2942adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
295c8b0795cSMark F. Adams    Input/Output Parameter:
2960cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
297c8b0795cSMark F. Adams */
298c8b0795cSMark F. Adams #undef __FUNCT__
299c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3002adfe9d3SBarry Smith static PetscErrorCode smoothAggs(Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
301c8b0795cSMark F. Adams {
302c8b0795cSMark F. Adams   PetscErrorCode ierr;
303c8b0795cSMark F. Adams   PetscBool      isMPI;
304c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
3053b4367a7SBarry Smith   MPI_Comm       comm;
3060cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
307c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
308c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3090cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3100cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
311c8b0795cSMark F. Adams   NState         *lid_state;
3120cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
313c8b0795cSMark F. Adams 
314c8b0795cSMark F. Adams   PetscFunctionBegin;
3153b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
316c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
317c8b0795cSMark F. Adams 
3180cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
319c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
320c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3213b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
322c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
323c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
324c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
325c8b0795cSMark F. Adams   }
326c8b0795cSMark F. Adams 
327c8b0795cSMark F. Adams   /* get submatrices */
328251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
329c8b0795cSMark F. Adams   if (isMPI) {
330c8b0795cSMark F. Adams     /* grab matrix objects */
331c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
332c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
333c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
334c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
335c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
336c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
337c8b0795cSMark F. Adams 
338c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
33911e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
340c8b0795cSMark F. Adams 
341785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
3420cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
343c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
344c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
345c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
346c8b0795cSMark F. Adams     }
347806fa848SBarry Smith   } else {
348c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
349c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3500298fd71SBarry Smith     lid_cprowID_1 = NULL;
351c8b0795cSMark F. Adams   }
35278a438d6SMark Adams   if (nloc>0) {
35371959b99SBarry Smith     if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
35471959b99SBarry Smith     if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
35571959b99SBarry Smith     if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
35671959b99SBarry Smith     if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
35778a438d6SMark Adams   }
358c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
359785e854fSJed Brown   ierr = PetscMalloc1(nloc, &lid_state);CHKERRQ(ierr);
360785e854fSJed Brown   ierr = PetscMalloc1(nloc, &lid_parent_gid);CHKERRQ(ierr);
361c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3620cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
363c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
364c8b0795cSMark F. Adams   }
3650cbbd2e1SMark F. Adams 
3660cbbd2e1SMark F. Adams   /* set lid_state */
3670cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
36841b27cdeSMark F. Adams     PetscCDPos pos;
369e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
370e78576d6SMark F. Adams     if (pos) {
371e78576d6SMark F. Adams       PetscInt gid1;
3722fa5cd67SKarl Rupp 
37371959b99SBarry Smith       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
37471959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3750cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
376b43b03e9SMark F. Adams     }
377b43b03e9SMark F. Adams   }
3780cbbd2e1SMark F. Adams 
3790cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
380c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
381c8b0795cSMark F. Adams     NState state = lid_state[lid];
382c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
38341b27cdeSMark F. Adams       PetscCDPos pos;
384e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
385e78576d6SMark F. Adams       while (pos) {
386e78576d6SMark F. Adams         PetscInt gid1;
387ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
388e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
389e78576d6SMark F. Adams 
3902fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
391c8b0795cSMark F. Adams       }
3920cbbd2e1SMark F. Adams     }
3930cbbd2e1SMark F. Adams   }
3940cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
395c8b0795cSMark F. Adams   if (isMPI) {
396c8b0795cSMark F. Adams     Vec tempVec;
397c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3982a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
399c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
400c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
401c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
402c8b0795cSMark F. Adams     }
403c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
404c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
405806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
406806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
407c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
408c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
409806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
411c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4120cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4130cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4140cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4150cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4160cbbd2e1SMark F. Adams     }
4170cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4180cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4190cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
420806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
421806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4220cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4230cbbd2e1SMark F. Adams 
424c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
425c8b0795cSMark F. Adams   } /* ismpi */
426c8b0795cSMark F. Adams 
427c8b0795cSMark F. Adams   /* doit */
428c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
429c8b0795cSMark F. Adams     NState state = lid_state[lid];
4300cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4310cbbd2e1SMark F. Adams       /* steal locals */
432c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
433c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
434c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4350cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
436c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4370cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4380cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4390cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4400cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
4410298fd71SBarry Smith             PetscCDPos pos,last=NULL;
442c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
443e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
444e78576d6SMark F. Adams             while (pos) {
445e78576d6SMark F. Adams               PetscInt gid;
446ffc955d6SMark F. Adams               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4470cbbd2e1SMark F. Adams               if (gid == gidj) {
44871959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
44941b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
45041b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4510cbbd2e1SMark F. Adams                 hav  = 1;
452c8b0795cSMark F. Adams                 break;
453806fa848SBarry Smith               } else last = pos;
454e78576d6SMark F. Adams 
455e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
456c8b0795cSMark F. Adams             }
457c8b0795cSMark F. Adams             if (hav!=1) {
45871959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
459c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
460c8b0795cSMark F. Adams             }
461806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4622fa5cd67SKarl Rupp             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 0.0' if the matrix is structurally symmetric.");
46341b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
464c8b0795cSMark F. Adams           }
465c8b0795cSMark F. Adams         }
4660cbbd2e1SMark F. Adams       } /* local neighbors */
467806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4680cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
469c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
470c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
471c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
472c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
473c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
474c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
475c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
476c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4770cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4780cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4790cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
4800298fd71SBarry Smith               PetscCDPos pos,last=NULL;
4810cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
482e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
483e78576d6SMark F. Adams               while (pos) {
484e78576d6SMark F. Adams                 PetscInt gid;
485ffc955d6SMark F. Adams                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4860cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4870cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
48871959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
48941b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4900cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4910cbbd2e1SMark F. Adams                   hav = 1;
492c8b0795cSMark F. Adams                   break;
493806fa848SBarry Smith                 } else last = pos;
494e78576d6SMark F. Adams 
495e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
496c8b0795cSMark F. Adams               }
497c8b0795cSMark F. Adams               if (hav!=1) {
498c666adbfSMark F. Adams                 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
499c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
500c8b0795cSMark F. Adams               }
501806fa848SBarry Smith             } else {
5020cbbd2e1SMark F. Adams               /* ghosts remove this later */
5030cbbd2e1SMark F. Adams             }
504c8b0795cSMark F. Adams           }
505c8b0795cSMark F. Adams         }
506c8b0795cSMark F. Adams       }
507c8b0795cSMark F. Adams     } /* selected/deleted */
508c8b0795cSMark F. Adams   } /* node loop */
509c8b0795cSMark F. Adams 
510c8b0795cSMark F. Adams   if (isMPI) {
5110cbbd2e1SMark F. Adams     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
5120cbbd2e1SMark F. Adams     Vec           tempVec,ghostgids2,ghostparents2;
5130cbbd2e1SMark F. Adams     PetscInt      cpid,nghost_2;
5140cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
515c8b0795cSMark F. Adams 
5160cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
5172a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5180cbbd2e1SMark F. Adams 
5190cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
520c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5210cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
522c8b0795cSMark F. Adams     }
523c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
524c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5250cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
526806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5280cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5290cbbd2e1SMark F. Adams 
5300cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5310cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5320cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5330cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5340cbbd2e1SMark F. Adams     }
5350cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5360cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5370cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
538806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5400cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5410cbbd2e1SMark F. Adams 
542c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
543c8b0795cSMark F. Adams 
5440cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
545f10ff945SMark F. Adams     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5460cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5470cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5480cbbd2e1SMark F. Adams       if (state==DELETED) {
5490cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5500cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5510cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5520cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5530cbbd2e1SMark F. Adams           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5540cbbd2e1SMark F. Adams         }
5550cbbd2e1SMark F. Adams       }
5560cbbd2e1SMark F. Adams     }
557c8b0795cSMark F. Adams 
5580cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
559c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
560c8b0795cSMark F. Adams       NState state = lid_state[lid];
561c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
5620298fd71SBarry Smith         PetscCDPos pos,last=NULL;
563c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
564e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
565e78576d6SMark F. Adams         while (pos) {
566e78576d6SMark F. Adams           PetscInt gid;
567ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
568e78576d6SMark F. Adams 
5690cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5700cbbd2e1SMark F. Adams             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5710cbbd2e1SMark F. Adams             if (cpid != -1) {
5720cbbd2e1SMark F. Adams               /* a moved ghost - */
5730cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
57441b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
575806fa848SBarry Smith             } else last = pos;
576806fa848SBarry Smith           } else last = pos;
577e78576d6SMark F. Adams 
578e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
579c8b0795cSMark F. Adams         } /* loop over list of deleted */
580c8b0795cSMark F. Adams       } /* selected */
581c8b0795cSMark F. Adams     }
5820cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
583c8b0795cSMark F. Adams 
5840cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5850cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5860cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5870cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5880cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5890cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
59041b27cdeSMark F. Adams         PetscCDPos pos;
5910cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
592e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
593e78576d6SMark F. Adams         while (pos) {
594e78576d6SMark F. Adams           PetscInt gidj;
595ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
596e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
597e78576d6SMark F. Adams 
5980cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
599c8b0795cSMark F. Adams         }
600c8b0795cSMark F. Adams         if (hav != 1) {
601ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
60241b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
603c8b0795cSMark F. Adams         }
604c8b0795cSMark F. Adams       }
605c8b0795cSMark F. Adams     }
606c8b0795cSMark F. Adams 
6070cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6080cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6090cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6100cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
611c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6120cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6130cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6140cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
615c8b0795cSMark F. Adams   }
616c8b0795cSMark F. Adams 
6170cbbd2e1SMark F. Adams   ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr);
618c8b0795cSMark F. Adams   ierr = PetscFree(lid_state);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"
634b8cd405aSMark F. Adams 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;
640b8cd405aSMark F. Adams 
6412e68589bSMark F. Adams   PetscFunctionBegin;
642b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
643b8cd405aSMark F. Adams   if (!mnull) {
644a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6459d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
64671959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
64771959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6480298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
649806fa848SBarry Smith   } else {
650b8cd405aSMark F. Adams     PetscReal *nullvec;
651b8cd405aSMark F. Adams     PetscBool has_const;
652b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
653b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6540298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
655b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
656a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
657785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
658b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
659b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
660b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
661b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
662b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
663b8cd405aSMark F. Adams     }
664b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
665b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6662fa5cd67SKarl Rupp 
66706e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
6682fa5cd67SKarl Rupp 
669b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
670b8cd405aSMark F. Adams   }
6712e68589bSMark F. Adams   PetscFunctionReturn(0);
6722e68589bSMark F. Adams }
6732e68589bSMark F. Adams 
6742e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6752e68589bSMark F. Adams /*
6762e68589bSMark F. Adams  formProl0
6772e68589bSMark F. Adams 
6782e68589bSMark F. Adams    Input Parameter:
6792adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
6802adfe9d3SBarry Smith    . bs - row block size
6810cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6820cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6832adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
6842e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6852e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6862e68589bSMark F. Adams   Output Parameter:
6872e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6882e68589bSMark F. Adams    . a_Prol - prolongation operator
6892e68589bSMark F. Adams */
6902e68589bSMark F. Adams #undef __FUNCT__
6912e68589bSMark F. Adams #define __FUNCT__ "formProl0"
6922adfe9d3SBarry 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)
6932e68589bSMark F. Adams {
6942e68589bSMark F. Adams   PetscErrorCode ierr;
6950cbbd2e1SMark F. Adams   PetscInt       Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
6963b4367a7SBarry Smith   MPI_Comm       comm;
69773911c69SBarry Smith   PetscMPIInt    rank;
6982e68589bSMark F. Adams   PetscReal      *out_data;
69941b27cdeSMark F. Adams   PetscCDPos     pos;
7000cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7010cbbd2e1SMark F. Adams 
702797e13b7SMark F. Adams /* #define OUT_AGGS */
703519f805aSKarl Rupp #if defined(OUT_AGGS)
7040298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7059057884aSMark F. Adams #endif
7062e68589bSMark F. Adams 
7072e68589bSMark F. Adams   PetscFunctionBegin;
7083b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7093b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7102e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
71171959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
71271959b99SBarry 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);
7130cbbd2e1SMark F. Adams   Iend   /= bs;
7140cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7152e68589bSMark F. Adams 
716f10ff945SMark F. Adams   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7170cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7180cbbd2e1SMark F. Adams     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7192e68589bSMark F. Adams   }
7202e68589bSMark F. Adams 
721519f805aSKarl Rupp #if defined(OUT_AGGS)
722c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7232fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7240cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7250cbbd2e1SMark F. Adams #endif
7260cbbd2e1SMark F. Adams 
7270cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7280cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
729e78576d6SMark F. Adams     PetscBool ise;
730e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
731e78576d6SMark F. Adams     if (!ise) nSelected++;
7320cbbd2e1SMark F. Adams   }
7330cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
73471959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
73571959b99SBarry 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);
7360cbbd2e1SMark F. Adams 
7372e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7380cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7392fa5cd67SKarl Rupp 
740785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
74133812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
7420cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7432e68589bSMark F. Adams 
7442e68589bSMark F. Adams   /* find points and set prolongation */
745c8b0795cSMark F. Adams   minsz = 100;
7462e68589bSMark F. Adams   ndone = 0;
7470cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
748e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
749e78576d6SMark F. Adams     if (jj > 0) {
7500cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7510cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7520cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7532e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7542e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7552e68589bSMark F. Adams       PetscInt       *fids;
75665d7b583SSatish Balay       PetscReal      *data;
7570cbbd2e1SMark F. Adams       /* count agg */
7580cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7590cbbd2e1SMark F. Adams 
7600cbbd2e1SMark F. Adams       /* get block */
761854ce69bSBarry Smith       ierr = PetscMalloc1(Mdata*N, &qqc);CHKERRQ(ierr);
762854ce69bSBarry Smith       ierr = PetscMalloc1(M*N, &qqr);CHKERRQ(ierr);
763785e854fSJed Brown       ierr = PetscMalloc1(N, &TAU);CHKERRQ(ierr);
764785e854fSJed Brown       ierr = PetscMalloc1(LWORK, &WORK);CHKERRQ(ierr);
765785e854fSJed Brown       ierr = PetscMalloc1(M, &fids);CHKERRQ(ierr);
7662e68589bSMark F. Adams 
7672e68589bSMark F. Adams       aggID = 0;
768e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
769e78576d6SMark F. Adams       while (pos) {
770e78576d6SMark F. Adams         PetscInt gid1;
771ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
772e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
773e78576d6SMark F. Adams 
7740cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7750cbbd2e1SMark F. Adams         else {
7760cbbd2e1SMark F. Adams           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
77771959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7780cbbd2e1SMark F. Adams         }
7792e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
78065d7b583SSatish Balay         data = &data_in[flid*bs];
7812e68589bSMark F. Adams         for (kk = ii = 0; ii < bs; ii++) {
7822e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7830cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7840cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7852e68589bSMark F. Adams           }
7862e68589bSMark F. Adams         }
787519f805aSKarl Rupp #if defined(OUT_AGGS)
788b2a4f308SMark F. Adams         if (llev==1) {
7899057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
7900cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
791c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
792f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
7930cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
794b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
795b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
7969057884aSMark F. Adams         }
7979057884aSMark F. Adams #endif
7982e68589bSMark F. Adams         /* set fine IDs */
7992e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8002e68589bSMark F. Adams 
8012e68589bSMark F. Adams         aggID++;
8020cbbd2e1SMark F. Adams       }
8032e68589bSMark F. Adams 
8042e68589bSMark F. Adams       /* pad with zeros */
8052e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8062e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8072e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8082e68589bSMark F. Adams         }
8092e68589bSMark F. Adams       }
8102e68589bSMark F. Adams 
8112e68589bSMark F. Adams       ndone += aggID;
8122e68589bSMark F. Adams       /* QR */
81384df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8148b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
81584df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
816d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8172e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8182e68589bSMark F. Adams       {
8192e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8202e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8212e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
82233812677SSatish 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);
8230cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8240cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8252e68589bSMark F. Adams           }
8262e68589bSMark F. Adams         }
8272e68589bSMark F. Adams       }
8282e68589bSMark F. Adams 
8292e68589bSMark F. Adams       /* get Q - row oriented */
8308b83055fSJed Brown       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
831c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8322e68589bSMark F. Adams 
8332e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8342e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8352e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8362e68589bSMark F. Adams         }
8372e68589bSMark F. Adams       }
8382e68589bSMark F. Adams 
8392e68589bSMark F. Adams       /* add diagonal block of P0 */
840c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
841c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
842c8b0795cSMark F. Adams       }
8432e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8442e68589bSMark F. Adams 
8452e68589bSMark F. Adams       ierr = PetscFree(qqc);CHKERRQ(ierr);
8462e68589bSMark F. Adams       ierr = PetscFree(qqr);CHKERRQ(ierr);
8472e68589bSMark F. Adams       ierr = PetscFree(TAU);CHKERRQ(ierr);
8482e68589bSMark F. Adams       ierr = PetscFree(WORK);CHKERRQ(ierr);
8492e68589bSMark F. Adams       ierr = PetscFree(fids);CHKERRQ(ierr);
850b43b03e9SMark F. Adams       clid++;
8510cbbd2e1SMark F. Adams     } /* coarse agg */
8520cbbd2e1SMark F. Adams   } /* for all fine nodes */
8530cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8540cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8552e68589bSMark F. Adams 
8563b4367a7SBarry Smith /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8572e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
8583b4367a7SBarry Smith /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8593b4367a7SBarry 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); */
8602e68589bSMark F. Adams 
861519f805aSKarl Rupp #if defined(OUT_AGGS)
862b2a4f308SMark F. Adams   if (llev==1) fclose(file);
8639057884aSMark F. Adams #endif
8640cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
8652e68589bSMark F. Adams   PetscFunctionReturn(0);
8662e68589bSMark F. Adams }
8672e68589bSMark F. Adams 
8685adeb434SBarry Smith #undef __FUNCT__
8695adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG_AGG"
8705adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
8715adeb434SBarry Smith {
8725adeb434SBarry Smith   PetscErrorCode ierr;
8735adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
8745adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
8755adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8765adeb434SBarry Smith 
8775adeb434SBarry Smith   PetscFunctionBegin;
8785adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
8795adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
8805adeb434SBarry Smith   PetscFunctionReturn(0);
8815adeb434SBarry Smith }
8825adeb434SBarry Smith 
8832e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
8842e68589bSMark F. Adams /*
885fd1112cbSBarry Smith    PCGAMGGraph_AGG
8862e68589bSMark F. Adams 
8872e68589bSMark F. Adams   Input Parameter:
8882e68589bSMark F. Adams    . pc - this
8892e68589bSMark F. Adams    . Amat - matrix on this fine level
8902e68589bSMark F. Adams   Output Parameter:
891c8b0795cSMark F. Adams    . a_Gmat -
8922e68589bSMark F. Adams */
8932e68589bSMark F. Adams #undef __FUNCT__
894fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_AGG"
8952adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
896c8b0795cSMark F. Adams {
897c8b0795cSMark F. Adams   PetscErrorCode            ierr;
898c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
899c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
900c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
901c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
902e0940f08SMark F. Adams   Mat                       Gmat;
9033b4367a7SBarry Smith   MPI_Comm                  comm;
90467744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
905c8b0795cSMark F. Adams 
906c8b0795cSMark F. Adams   PetscFunctionBegin;
9073b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
908fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
909c8b0795cSMark F. Adams 
91067744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
911c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9120cbbd2e1SMark F. Adams 
9132d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
914bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
915c8b0795cSMark F. Adams 
916e0940f08SMark F. Adams   *a_Gmat = Gmat;
917c8b0795cSMark F. Adams 
918fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
919c8b0795cSMark F. Adams   PetscFunctionReturn(0);
920c8b0795cSMark F. Adams }
921c8b0795cSMark F. Adams 
922c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
923c8b0795cSMark F. Adams /*
924b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
925c8b0795cSMark F. Adams 
926c8b0795cSMark F. Adams   Input Parameter:
927e0940f08SMark F. Adams    . a_pc - this
928e0940f08SMark F. Adams   Input/Output Parameter:
9290cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
930c8b0795cSMark F. Adams   Output Parameter:
9310cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
932c8b0795cSMark F. Adams */
933c8b0795cSMark F. Adams #undef __FUNCT__
934b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
9352fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
936c8b0795cSMark F. Adams {
937c8b0795cSMark F. Adams   PetscErrorCode ierr;
938e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
939c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
940c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9410cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9420cbbd2e1SMark F. Adams   IS             perm;
943c8b0795cSMark F. Adams   PetscInt       Ii,nloc,bs,n,m;
944c8b0795cSMark F. Adams   PetscInt       *permute;
945c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
946b43b03e9SMark F. Adams   MatCoarsen     crs;
9473b4367a7SBarry Smith   MPI_Comm       comm;
94873911c69SBarry Smith   PetscMPIInt    rank;
949c8b0795cSMark F. Adams 
950c8b0795cSMark F. Adams   PetscFunctionBegin;
9510cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9523b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9533b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
954e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
95571959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
95671959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
957c8b0795cSMark F. Adams   nloc = n/bs;
958c8b0795cSMark F. Adams 
95957d29afaSToby Isaac   if (pc_gamg->firstCoarsen && pc_gamg_agg->square_graph) { /* we don't want to square coarse grids */
960878e152fSMark F. Adams     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
961806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
962806fa848SBarry Smith   } else Gmat2 = Gmat1;
963c8b0795cSMark F. Adams 
964c8b0795cSMark F. Adams   /* get MIS aggs */
965c8b0795cSMark F. Adams   /* randomize */
966785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
967785e854fSJed Brown   ierr = PetscMalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
968c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
969c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
970c8b0795cSMark F. Adams     permute[Ii]   = Ii;
971c8b0795cSMark F. Adams   }
972c8b0795cSMark F. Adams   srand(1); /* make deterministic */
973c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
974c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
975c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
976c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
977c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
978c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
979c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
980c8b0795cSMark F. Adams     }
981c8b0795cSMark F. Adams   }
982c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
983c8b0795cSMark F. Adams 
984806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
9850cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9860cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
987b43b03e9SMark F. Adams #endif
9883b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
9899057884aSMark F. Adams   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
9909057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
991b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
992b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
993b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
994b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
9950cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
996b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
997b43b03e9SMark F. Adams 
998c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
999c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10000cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10010cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1002b43b03e9SMark F. Adams #endif
10039f3f12c8SMark F. Adams 
1004c8b0795cSMark F. Adams   /* smooth aggs */
1005e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10060cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10070cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1008c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1009e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
101041b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10113b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1012806fa848SBarry Smith   } else {
10130cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1014f10ff945SMark F. Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
101541b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10160cbbd2e1SMark F. Adams     if (mat) {
10170cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10180cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10190cbbd2e1SMark F. Adams     }
10200cbbd2e1SMark F. Adams   }
10210cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1022c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1023c8b0795cSMark F. Adams }
1024c8b0795cSMark F. Adams 
1025c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1026c8b0795cSMark F. Adams /*
10270cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1028c8b0795cSMark F. Adams 
1029c8b0795cSMark F. Adams  Input Parameter:
1030c8b0795cSMark F. Adams  . pc - this
1031c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1032c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10330cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1034c8b0795cSMark F. Adams  Output Parameter:
1035c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1036c8b0795cSMark F. Adams  */
1037c8b0795cSMark F. Adams #undef __FUNCT__
10380cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
10392adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10402e68589bSMark F. Adams {
10412e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10422e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
1043c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10442e68589bSMark F. Adams   PetscErrorCode ierr;
1045c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1046c8b0795cSMark F. Adams   Mat            Prol;
1047c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10483b4367a7SBarry Smith   MPI_Comm       comm;
10490cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1050c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1051c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1052*d9558ea9SBarry Smith   MatType        mtype;
10532e68589bSMark F. Adams 
10542e68589bSMark F. Adams   PetscFunctionBegin;
10553b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
10560cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10573b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10583b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10592e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1060c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
106171959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
106271959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
10632e68589bSMark F. Adams 
10642e68589bSMark F. Adams   /* get 'nLocalSelected' */
10650cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1066e78576d6SMark F. Adams     PetscBool ise;
10670cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1068e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1069e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
10702e68589bSMark F. Adams   }
10712e68589bSMark F. Adams 
10722e68589bSMark F. Adams   /* create prolongator, create P matrix */
1073*d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
10743b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1075806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1076a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1077*d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
10780298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
10790298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1080a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1081a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
10820298fd71SBarry Smith   /* data_cols, NULL, data_cols, NULL, */
1083a2f3521dSMark F. Adams   /* &Prol); */
10842e68589bSMark F. Adams 
10852e68589bSMark F. Adams   /* can get all points "removed" */
1086c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1087c8b0795cSMark F. Adams   if (ii==0) {
1088bf4339c2SBarry Smith     ierr = PetscInfo(pc,"no selected points on coarse grid\n");CHKERRQ(ierr);
10892e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
10900298fd71SBarry Smith     *a_P_out = NULL;  /* out */
10912e68589bSMark F. Adams     PetscFunctionReturn(0);
10922e68589bSMark F. Adams   }
1093bf4339c2SBarry Smith   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1094c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
10950cbbd2e1SMark F. Adams 
109671959b99SBarry 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);
1097c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
109871959b99SBarry 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);
10992e68589bSMark F. Adams 
11002e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11010cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11020cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11032e68589bSMark F. Adams #endif
1104c5df96a5SBarry Smith   if (size > 1) { /*  */
11052e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1106785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
11072e68589bSMark F. Adams     for (jj = 0; jj < data_cols; jj++) {
1108c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1109a2f3521dSMark F. Adams         PetscInt        ii,stride;
1110c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11112fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11122fa5cd67SKarl Rupp 
1113806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1114a2f3521dSMark F. Adams 
11152e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1116785e854fSJed Brown           ierr    = PetscMalloc1(stride*bs*data_cols, &data_w_ghost);CHKERRQ(ierr);
1117a2f3521dSMark F. Adams           nbnodes = bs*stride;
11182e68589bSMark F. Adams         }
1119a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11202fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11212e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11222e68589bSMark F. Adams       }
11232e68589bSMark F. Adams     }
11242e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1125806fa848SBarry Smith   } else {
1126c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1127c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11282e68589bSMark F. Adams   }
11292e68589bSMark F. Adams 
11302e68589bSMark F. Adams   /* get P0 */
1131c5df96a5SBarry Smith   if (size > 1) {
11322e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1133a2f3521dSMark F. Adams     PetscInt  stride;
11342e68589bSMark F. Adams 
1135785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
11362e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1137806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1138785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1139a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11402e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1141a2f3521dSMark F. Adams 
114271959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11432e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1144806fa848SBarry Smith   } else {
1145785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
11462e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11472e68589bSMark F. Adams   }
11480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11490cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11500cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
11512e68589bSMark F. Adams #endif
1152c8b0795cSMark F. Adams   {
11530298fd71SBarry Smith     PetscReal *data_out = NULL;
1154806fa848SBarry Smith     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1155c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11562fa5cd67SKarl Rupp 
1157c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
1158c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1159c8b0795cSMark F. Adams     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1160c8b0795cSMark F. Adams   }
11610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11620cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1163c8b0795cSMark F. Adams #endif
1164c5df96a5SBarry Smith   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
11652e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
11662e68589bSMark F. Adams 
1167c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1168ed4e82eaSPeter Brune 
11690cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1170c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1171c8b0795cSMark F. Adams }
1172c8b0795cSMark F. Adams 
1173c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1174c8b0795cSMark F. Adams /*
1175fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1176c8b0795cSMark F. Adams 
1177c8b0795cSMark F. Adams   Input Parameter:
1178c8b0795cSMark F. Adams    . pc - this
1179c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1180c8b0795cSMark F. Adams  In/Output Parameter:
1181c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1182c8b0795cSMark F. Adams */
1183c8b0795cSMark F. Adams #undef __FUNCT__
1184fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_AGG"
11852adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1186c8b0795cSMark F. Adams {
1187c8b0795cSMark F. Adams   PetscErrorCode ierr;
1188c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1189c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1190c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1191c8b0795cSMark F. Adams   PetscInt       jj;
1192c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
11933b4367a7SBarry Smith   MPI_Comm       comm;
1194c8b0795cSMark F. Adams 
1195c8b0795cSMark F. Adams   PetscFunctionBegin;
11963b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1197fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1198c8b0795cSMark F. Adams 
11992e68589bSMark F. Adams   /* smooth P0 */
1200c8b0795cSMark F. Adams   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12012e68589bSMark F. Adams     Mat       tMat;
12022e68589bSMark F. Adams     Vec       diag;
12032e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12040cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12050cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12062e68589bSMark F. Adams #endif
12072e68589bSMark F. Adams     if (jj == 0) {
12082e68589bSMark F. Adams       KSP eksp;
12092e68589bSMark F. Adams       Vec bb, xx;
1210da509fc8SJed Brown       PC  epc;
12112a7a6963SBarry Smith       ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
12122a7a6963SBarry Smith       ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
12132e68589bSMark F. Adams       {
12142e68589bSMark F. Adams         PetscRandom rctx;
12153b4367a7SBarry Smith         ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
12162e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12172e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12182e68589bSMark F. Adams         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
12192e68589bSMark F. Adams       }
1220e696ed0bSMark F. Adams 
1221e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1222e696ed0bSMark F. Adams       {
1223e696ed0bSMark F. Adams         PetscInt    Istart,Iend,ncols,jj,Ii;
1224e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1225e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1226e696ed0bSMark F. Adams         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1227e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1228e696ed0bSMark F. Adams           if (ncols <= 1) {
1229e696ed0bSMark F. Adams             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1230e696ed0bSMark F. Adams           }
1231e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1232e696ed0bSMark F. Adams         }
1233e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1234e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1235e696ed0bSMark F. Adams       }
1236e696ed0bSMark F. Adams 
12373b4367a7SBarry Smith       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1238806fa848SBarry Smith       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12392e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1240da509fc8SJed Brown       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1241c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1242c2436cacSMark F. Adams       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1243c2436cacSMark F. Adams 
1244c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
124523ee1639SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
12462e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
12472e68589bSMark F. Adams 
1248da509fc8SJed Brown       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1249da509fc8SJed Brown       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1250c2436cacSMark F. Adams 
12512e68589bSMark F. Adams       /* solve - keep stuff out of logging */
12522e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
12532e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
12542e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
12552e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
12562e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
12572e68589bSMark F. Adams 
12582e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1259bf4339c2SBarry Smith       ierr = PetscInfo3(pc,"\t\t\t smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
12602e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
12612e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
12622e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
12632e68589bSMark F. Adams 
12642e68589bSMark F. Adams       if (pc_gamg->emax_id == -1) {
12653bf036e2SBarry Smith         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
126671959b99SBarry Smith         if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1");
12672e68589bSMark F. Adams       }
12682e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
12692e68589bSMark F. Adams     }
12702e68589bSMark F. Adams 
12712e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
12722e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
12732a7a6963SBarry Smith     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
12742e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
12752e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
12762e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
12772e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1278b4ec6429SMark F. Adams     alpha = -1.4/emax;
12792e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
12802e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
12812e68589bSMark F. Adams     Prol  = tMat;
12820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12830cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12842e68589bSMark F. Adams #endif
12852e68589bSMark F. Adams   }
1286fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1287c8b0795cSMark F. Adams   *a_P = Prol;
1288c8b0795cSMark F. Adams   PetscFunctionReturn(0);
12892e68589bSMark F. Adams }
12902e68589bSMark F. Adams 
1291c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1292c8b0795cSMark F. Adams /*
1293c8b0795cSMark F. Adams    PCCreateGAMG_AGG
12942e68589bSMark F. Adams 
1295c8b0795cSMark F. Adams   Input Parameter:
1296c8b0795cSMark F. Adams    . pc -
1297c8b0795cSMark F. Adams */
1298c8b0795cSMark F. Adams #undef __FUNCT__
1299c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1300c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1301c8b0795cSMark F. Adams {
1302c8b0795cSMark F. Adams   PetscErrorCode ierr;
1303c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1304c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1305c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
13062e68589bSMark F. Adams 
1307c8b0795cSMark F. Adams   PetscFunctionBegin;
1308c8b0795cSMark F. Adams   /* create sub context for SA */
1309b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1310c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1311c8b0795cSMark F. Adams 
13121ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
13139b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1314c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1315c8b0795cSMark F. Adams 
1316c8b0795cSMark F. Adams   /* set internal function pointers */
1317fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
13181ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
13191ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1320fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
13211ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
13225adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1323c8b0795cSMark F. Adams 
1324bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
13252e68589bSMark F. Adams   PetscFunctionReturn(0);
13262e68589bSMark F. Adams }
1327