xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 1ab5ffc991220836b43ae139cc95ce5cf8f2cb1d)
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"
1512e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc)
1522e68589bSMark F. Adams {
1532e68589bSMark F. Adams   PetscErrorCode ierr;
1542e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1552e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
156c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1572e68589bSMark F. Adams   PetscBool      flag;
1582e68589bSMark F. Adams 
1592e68589bSMark F. Adams   PetscFunctionBegin;
1602e68589bSMark F. Adams 
1612e68589bSMark F. Adams   ierr = PetscOptionsHead("GAMG-AGG options");CHKERRQ(ierr);
1622e68589bSMark F. Adams   {
1632e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
164c8b0795cSMark F. Adams     pc_gamg_agg->nsmooths = 0;
1652fa5cd67SKarl Rupp 
1662e68589bSMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
1672e68589bSMark F. Adams                            "smoothing steps for smoothed aggregation, usually 1 (0)",
1682e68589bSMark F. Adams                            "PCGAMGSetNSmooths",
169c8b0795cSMark F. Adams                            pc_gamg_agg->nsmooths,
170c8b0795cSMark F. Adams                            &pc_gamg_agg->nsmooths,
171806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
172c8b0795cSMark F. Adams 
173c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
174c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
1752fa5cd67SKarl Rupp 
176c8b0795cSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_sym_graph",
177581a99e3SJed Brown                             "Set for asymmetric matrices",
178c8b0795cSMark F. Adams                             "PCGAMGSetSymGraph",
179c8b0795cSMark F. Adams                             pc_gamg_agg->sym_graph,
180c8b0795cSMark F. Adams                             &pc_gamg_agg->sym_graph,
181806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
182ef4ad70eSMark F. Adams 
183ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
184ef4ad70eSMark F. Adams     pc_gamg_agg->square_graph = PETSC_TRUE;
1852fa5cd67SKarl Rupp 
186ef4ad70eSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_square_graph",
1870cbbd2e1SMark F. Adams                             "For faster coarsening and lower coarse grid complexity",
188ef4ad70eSMark F. Adams                             "PCGAMGSetSquareGraph",
189ef4ad70eSMark F. Adams                             pc_gamg_agg->square_graph,
190ef4ad70eSMark F. Adams                             &pc_gamg_agg->square_graph,
191806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
1922e68589bSMark F. Adams   }
1932e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1942e68589bSMark F. Adams   PetscFunctionReturn(0);
1952e68589bSMark F. Adams }
1962e68589bSMark F. Adams 
1972e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1982e68589bSMark F. Adams /*
1992e68589bSMark F. Adams    PCDestroy_AGG
2002e68589bSMark F. Adams 
2012e68589bSMark F. Adams   Input Parameter:
2022e68589bSMark F. Adams    . pc -
2032e68589bSMark F. Adams */
2042e68589bSMark F. Adams #undef __FUNCT__
2052e68589bSMark F. Adams #define __FUNCT__ "PCDestroy_AGG"
2062e68589bSMark F. Adams PetscErrorCode PCDestroy_AGG(PC pc)
2072e68589bSMark F. Adams {
2082e68589bSMark F. Adams   PetscErrorCode ierr;
2092e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
2102e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
211c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
2122e68589bSMark F. Adams 
2132e68589bSMark F. Adams   PetscFunctionBegin;
214c8b0795cSMark F. Adams   if (pc_gamg_agg) {
215c8b0795cSMark F. Adams     ierr        = PetscFree(pc_gamg_agg);CHKERRQ(ierr);
216c8b0795cSMark F. Adams     pc_gamg_agg = 0;
2172e68589bSMark F. Adams   }
2182e68589bSMark F. Adams 
2192e68589bSMark F. Adams   /* call base class */
2202e68589bSMark F. Adams   ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
2212e68589bSMark F. Adams   PetscFunctionReturn(0);
2222e68589bSMark F. Adams }
2232e68589bSMark F. Adams 
2242e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2252e68589bSMark F. Adams /*
2262e68589bSMark F. Adams    PCSetCoordinates_AGG
227302f38e8SMark F. Adams      - collective
2282e68589bSMark F. Adams 
2292e68589bSMark F. Adams    Input Parameter:
2302e68589bSMark F. Adams    . pc - the preconditioner context
231a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
232302f38e8SMark F. Adams    . a_nloc - number of vertices local
233302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2342e68589bSMark F. Adams */
2351e6b0712SBarry Smith 
2362e68589bSMark F. Adams #undef __FUNCT__
2372e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2381e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2392e68589bSMark F. Adams {
2402e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2412e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2422e68589bSMark F. Adams   PetscErrorCode ierr;
24369344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
244a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2452e68589bSMark F. Adams 
2462e68589bSMark F. Adams   PetscFunctionBegin;
247a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
248a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
249302f38e8SMark F. Adams   nloc = a_nloc;
2502e68589bSMark F. Adams 
2512e68589bSMark F. Adams   /* SA: null space vectors */
25269344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
25369344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
254a2f3521dSMark F. Adams   else if (coords) {
255c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
25669344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
25769344418SMark F. Adams     if (ndm != ndf) {
258c666adbfSMark 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);
259a2f3521dSMark F. Adams     }
260806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
26171959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
26271959b99SBarry 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);
263c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2642e68589bSMark F. Adams 
2652e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2662e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2672e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
268302f38e8SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
2692e68589bSMark F. Adams   }
2702e68589bSMark F. Adams   /* copy data in - column oriented */
2712e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
27269344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
27369344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
274c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2752e68589bSMark F. Adams     else {
27669344418SMark F. Adams       /* translational modes */
2772fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2782fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
27969344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2802e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2812fa5cd67SKarl Rupp         }
2822fa5cd67SKarl Rupp       }
28369344418SMark F. Adams 
28469344418SMark F. Adams       /* rotational modes */
2852e68589bSMark F. Adams       if (coords) {
28669344418SMark F. Adams         if (ndm == 2) {
2872e68589bSMark F. Adams           data   += 2*M;
2882e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2892e68589bSMark F. Adams           data[1] =  coords[2*kk];
290806fa848SBarry Smith         } else {
2912e68589bSMark F. Adams           data   += 3*M;
2922e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2932e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2942e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2952e68589bSMark F. Adams         }
2962e68589bSMark F. Adams       }
2972e68589bSMark F. Adams     }
2982e68589bSMark F. Adams   }
2992e68589bSMark F. Adams 
3002e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
3012e68589bSMark F. Adams   PetscFunctionReturn(0);
3022e68589bSMark F. Adams }
3032e68589bSMark F. Adams 
304b43b03e9SMark F. Adams typedef PetscInt NState;
305b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
306b43b03e9SMark F. Adams static const NState DELETED =-1;
307b43b03e9SMark F. Adams static const NState REMOVED =-3;
308b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
309b43b03e9SMark F. Adams 
310c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
311c8b0795cSMark F. Adams /*
312b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
313b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
314c8b0795cSMark F. Adams 
315c8b0795cSMark F. Adams    Input Parameter:
316c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
317c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
318c8b0795cSMark F. Adams    Input/Output Parameter:
3190cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
320c8b0795cSMark F. Adams */
321c8b0795cSMark F. Adams #undef __FUNCT__
322c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3230cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
3240cbbd2e1SMark F. Adams                                  const Mat Gmat_1,  /* base graph */
3250cbbd2e1SMark F. Adams                                  /* const IS selected_2, [nselected local] selected vertices */
3261147fc2aSKarl Rupp                                  PetscCoarsenData *aggs_2)  /* [nselected local] global ID of aggregate */
327c8b0795cSMark F. Adams {
328c8b0795cSMark F. Adams   PetscErrorCode ierr;
329c8b0795cSMark F. Adams   PetscBool      isMPI;
330c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
3313b4367a7SBarry Smith   MPI_Comm       comm;
332c5df96a5SBarry Smith   PetscMPIInt    rank,size;
3330cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
334c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
335c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3360cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3370cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
338c8b0795cSMark F. Adams   NState         *lid_state;
3390cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
340c8b0795cSMark F. Adams 
341c8b0795cSMark F. Adams   PetscFunctionBegin;
3423b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
3433b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
3443b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
345c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
346c8b0795cSMark F. Adams 
3470cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
348c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
349c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3503b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
351c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
352c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
353c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
354c8b0795cSMark F. Adams   }
355c8b0795cSMark F. Adams 
356c8b0795cSMark F. Adams   /* get submatrices */
357251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
358c8b0795cSMark F. Adams   if (isMPI) {
359c8b0795cSMark F. Adams     /* grab matrix objects */
360c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
361c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
362c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
363c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
364c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
365c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
366c8b0795cSMark F. Adams 
367c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
368c8b0795cSMark F. Adams     matB_1->compressedrow.check = PETSC_TRUE;
3692fa5cd67SKarl Rupp 
370806fa848SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
371c8b0795cSMark F. Adams 
372c8b0795cSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &lid_cprowID_1);CHKERRQ(ierr);
3730cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
374c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
375c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
376c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
377c8b0795cSMark F. Adams     }
378806fa848SBarry Smith   } else {
379c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
380c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3810298fd71SBarry Smith     lid_cprowID_1 = NULL;
382c8b0795cSMark F. Adams   }
38371959b99SBarry Smith   if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
38471959b99SBarry Smith   if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
38571959b99SBarry Smith   if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
38671959b99SBarry Smith   if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
387c8b0795cSMark F. Adams 
388c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
389c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(NState), &lid_state);CHKERRQ(ierr);
3900cbbd2e1SMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscScalar), &lid_parent_gid);CHKERRQ(ierr);
391c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3920cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
393c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
394c8b0795cSMark F. Adams   }
3950cbbd2e1SMark F. Adams 
3960cbbd2e1SMark F. Adams   /* set lid_state */
3970cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
39841b27cdeSMark F. Adams     PetscCDPos pos;
399e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
400e78576d6SMark F. Adams     if (pos) {
401e78576d6SMark F. Adams       PetscInt gid1;
4022fa5cd67SKarl Rupp 
40371959b99SBarry Smith       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
40471959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
4050cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
406b43b03e9SMark F. Adams     }
407b43b03e9SMark F. Adams   }
4080cbbd2e1SMark F. Adams 
4090cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
410c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
411c8b0795cSMark F. Adams     NState state = lid_state[lid];
412c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
41341b27cdeSMark F. Adams       PetscCDPos pos;
414e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
415e78576d6SMark F. Adams       while (pos) {
416e78576d6SMark F. Adams         PetscInt gid1;
417ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
418e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
419e78576d6SMark F. Adams 
4202fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
421c8b0795cSMark F. Adams       }
4220cbbd2e1SMark F. Adams     }
4230cbbd2e1SMark F. Adams   }
4240cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
425c8b0795cSMark F. Adams   if (isMPI) {
426c8b0795cSMark F. Adams     Vec tempVec;
427c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
428c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
429c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
430c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
431c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
432c8b0795cSMark F. Adams     }
433c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
434c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
435806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
436806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
437c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
438c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
439806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
440806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
441c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4420cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4430cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4440cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4450cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4460cbbd2e1SMark F. Adams     }
4470cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4480cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4490cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
450806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
451806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4520cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4530cbbd2e1SMark F. Adams 
454c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
455c8b0795cSMark F. Adams   } /* ismpi */
456c8b0795cSMark F. Adams 
457c8b0795cSMark F. Adams   /* doit */
458c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
459c8b0795cSMark F. Adams     NState state = lid_state[lid];
4600cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4610cbbd2e1SMark F. Adams       /* steal locals */
462c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
463c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
464c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4650cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
466c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4670cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4680cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4690cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4700cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
4710298fd71SBarry Smith             PetscCDPos pos,last=NULL;
472c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
473e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
474e78576d6SMark F. Adams             while (pos) {
475e78576d6SMark F. Adams               PetscInt gid;
476ffc955d6SMark F. Adams               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4770cbbd2e1SMark F. Adams               if (gid == gidj) {
47871959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
47941b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
48041b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4810cbbd2e1SMark F. Adams                 hav  = 1;
482c8b0795cSMark F. Adams                 break;
483806fa848SBarry Smith               } else last = pos;
484e78576d6SMark F. Adams 
485e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
486c8b0795cSMark F. Adams             }
487c8b0795cSMark F. Adams             if (hav!=1) {
48871959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
489c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
490c8b0795cSMark F. Adams             }
491806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4922fa5cd67SKarl 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.");
49341b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
494c8b0795cSMark F. Adams           }
495c8b0795cSMark F. Adams         }
4960cbbd2e1SMark F. Adams       } /* local neighbors */
497806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4980cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
499c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
500c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
501c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
502c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
503c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
504c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
505c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
506c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
5070cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
5080cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
5090cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
5100298fd71SBarry Smith               PetscCDPos pos,last=NULL;
5110cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
512e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
513e78576d6SMark F. Adams               while (pos) {
514e78576d6SMark F. Adams                 PetscInt gid;
515ffc955d6SMark F. Adams                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
5160cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
5170cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
51871959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
51941b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
5200cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
5210cbbd2e1SMark F. Adams                   hav = 1;
522c8b0795cSMark F. Adams                   break;
523806fa848SBarry Smith                 } else last = pos;
524e78576d6SMark F. Adams 
525e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
526c8b0795cSMark F. Adams               }
527c8b0795cSMark F. Adams               if (hav!=1) {
528c666adbfSMark F. Adams                 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
529c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
530c8b0795cSMark F. Adams               }
531806fa848SBarry Smith             } else {
5320cbbd2e1SMark F. Adams               /* ghosts remove this later */
5330cbbd2e1SMark F. Adams             }
534c8b0795cSMark F. Adams           }
535c8b0795cSMark F. Adams         }
536c8b0795cSMark F. Adams       }
537c8b0795cSMark F. Adams     } /* selected/deleted */
538c8b0795cSMark F. Adams   } /* node loop */
539c8b0795cSMark F. Adams 
540c8b0795cSMark F. Adams   if (isMPI) {
5410cbbd2e1SMark F. Adams     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
5420cbbd2e1SMark F. Adams     Vec           tempVec,ghostgids2,ghostparents2;
5430cbbd2e1SMark F. Adams     PetscInt      cpid,nghost_2;
5440cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
545c8b0795cSMark F. Adams 
5460cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
547c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5480cbbd2e1SMark F. Adams 
5490cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
550c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5510cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
552c8b0795cSMark F. Adams     }
553c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
554c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5550cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
556806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
557806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5580cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5590cbbd2e1SMark F. Adams 
5600cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5610cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5620cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5630cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5640cbbd2e1SMark F. Adams     }
5650cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5660cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5670cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
568806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
569806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5700cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5710cbbd2e1SMark F. Adams 
572c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
573c8b0795cSMark F. Adams 
5740cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
575f10ff945SMark F. Adams     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5760cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5770cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5780cbbd2e1SMark F. Adams       if (state==DELETED) {
5790cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5800cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5810cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5820cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5830cbbd2e1SMark F. Adams           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5840cbbd2e1SMark F. Adams         }
5850cbbd2e1SMark F. Adams       }
5860cbbd2e1SMark F. Adams     }
587c8b0795cSMark F. Adams 
5880cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
589c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
590c8b0795cSMark F. Adams       NState state = lid_state[lid];
591c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
5920298fd71SBarry Smith         PetscCDPos pos,last=NULL;
593c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
594e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
595e78576d6SMark F. Adams         while (pos) {
596e78576d6SMark F. Adams           PetscInt gid;
597ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
598e78576d6SMark F. Adams 
5990cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
6000cbbd2e1SMark F. Adams             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
6010cbbd2e1SMark F. Adams             if (cpid != -1) {
6020cbbd2e1SMark F. Adams               /* a moved ghost - */
6030cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
60441b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
605806fa848SBarry Smith             } else last = pos;
606806fa848SBarry Smith           } else last = pos;
607e78576d6SMark F. Adams 
608e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
609c8b0795cSMark F. Adams         } /* loop over list of deleted */
610c8b0795cSMark F. Adams       } /* selected */
611c8b0795cSMark F. Adams     }
6120cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
613c8b0795cSMark F. Adams 
6140cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
6150cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
6160cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6170cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
6180cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6190cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
62041b27cdeSMark F. Adams         PetscCDPos pos;
6210cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
622e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
623e78576d6SMark F. Adams         while (pos) {
624e78576d6SMark F. Adams           PetscInt gidj;
625ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
626e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
627e78576d6SMark F. Adams 
6280cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
629c8b0795cSMark F. Adams         }
630c8b0795cSMark F. Adams         if (hav != 1) {
631ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
63241b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
633c8b0795cSMark F. Adams         }
634c8b0795cSMark F. Adams       }
635c8b0795cSMark F. Adams     }
636c8b0795cSMark F. Adams 
6370cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6380cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6390cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6400cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
641c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6420cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6430cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6440cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
645c8b0795cSMark F. Adams   }
646c8b0795cSMark F. Adams 
6470cbbd2e1SMark F. Adams   ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr);
648c8b0795cSMark F. Adams   ierr = PetscFree(lid_state);CHKERRQ(ierr);
649c8b0795cSMark F. Adams   PetscFunctionReturn(0);
650c8b0795cSMark F. Adams }
6512e68589bSMark F. Adams 
6522e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6532e68589bSMark F. Adams /*
654a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
655a2f3521dSMark F. Adams       Looks in Mat for near null space.
656a2f3521dSMark F. Adams       Does not work for Stokes
6572e68589bSMark F. Adams 
6582e68589bSMark F. Adams   Input Parameter:
6592e68589bSMark F. Adams    . pc -
660a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6612e68589bSMark F. Adams */
6622e68589bSMark F. Adams #undef __FUNCT__
6632e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
664b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6652e68589bSMark F. Adams {
6662e68589bSMark F. Adams   PetscErrorCode ierr;
667b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
668b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
669b8cd405aSMark F. Adams   MatNullSpace   mnull;
670b8cd405aSMark F. Adams 
6712e68589bSMark F. Adams   PetscFunctionBegin;
672b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
673b8cd405aSMark F. Adams   if (!mnull) {
674a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6759d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
67671959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
67771959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6780298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
679806fa848SBarry Smith   } else {
680b8cd405aSMark F. Adams     PetscReal *nullvec;
681b8cd405aSMark F. Adams     PetscBool has_const;
682b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
683b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6840298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
685b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
6868caf3d72SBarry Smith     ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
687b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
688b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
689b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
690b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
691b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
692b8cd405aSMark F. Adams     }
693b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
694b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6952fa5cd67SKarl Rupp 
696a2f3521dSMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); /* this does not work for Stokes */
6972fa5cd67SKarl Rupp 
698b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
699b8cd405aSMark F. Adams   }
7002e68589bSMark F. Adams   PetscFunctionReturn(0);
7012e68589bSMark F. Adams }
7022e68589bSMark F. Adams 
7032e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
7042e68589bSMark F. Adams /*
7052e68589bSMark F. Adams  formProl0
7062e68589bSMark F. Adams 
7072e68589bSMark F. Adams    Input Parameter:
7080cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
7092e68589bSMark F. Adams    . bs - block size
7100cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7110cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7122e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
7132e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7142e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7152e68589bSMark F. Adams   Output Parameter:
7162e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7172e68589bSMark F. Adams    . a_Prol - prolongation operator
7182e68589bSMark F. Adams */
7192e68589bSMark F. Adams #undef __FUNCT__
7202e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7210cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */
7220cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
7230cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
7240cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
7250cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7260cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7270cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7282e68589bSMark F. Adams                                 PetscReal **a_data_out,
7291147fc2aSKarl Rupp                                 Mat a_Prol) /* prolongation operator (output)*/
7302e68589bSMark F. Adams {
7312e68589bSMark F. Adams   PetscErrorCode ierr;
7320cbbd2e1SMark F. Adams   PetscInt       Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7333b4367a7SBarry Smith   MPI_Comm       comm;
734c5df96a5SBarry Smith   PetscMPIInt    rank, size;
7352e68589bSMark F. Adams   PetscReal      *out_data;
73641b27cdeSMark F. Adams   PetscCDPos     pos;
7370cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7380cbbd2e1SMark F. Adams 
739797e13b7SMark F. Adams /* #define OUT_AGGS */
740519f805aSKarl Rupp #if defined(OUT_AGGS)
7410298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7429057884aSMark F. Adams #endif
7432e68589bSMark F. Adams 
7442e68589bSMark F. Adams   PetscFunctionBegin;
7453b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7463b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7473b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
7482e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
74971959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
75071959b99SBarry 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);
7510cbbd2e1SMark F. Adams   Iend   /= bs;
7520cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7532e68589bSMark F. Adams 
754f10ff945SMark F. Adams   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7550cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7560cbbd2e1SMark F. Adams     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7572e68589bSMark F. Adams   }
7582e68589bSMark F. Adams 
759519f805aSKarl Rupp #if defined(OUT_AGGS)
760c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7612fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7620cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7630cbbd2e1SMark F. Adams #endif
7640cbbd2e1SMark F. Adams 
7650cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7660cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
767e78576d6SMark F. Adams     PetscBool ise;
768e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
769e78576d6SMark F. Adams     if (!ise) nSelected++;
7700cbbd2e1SMark F. Adams   }
7710cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
77271959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
77371959b99SBarry 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);
7740cbbd2e1SMark F. Adams 
7752e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7760cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7772fa5cd67SKarl Rupp 
7780cbbd2e1SMark F. Adams   ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr);
7792fa5cd67SKarl Rupp   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=1.e300;
7800cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7812e68589bSMark F. Adams 
7822e68589bSMark F. Adams   /* find points and set prolongation */
783c8b0795cSMark F. Adams   minsz = 100;
7842e68589bSMark F. Adams   ndone = 0;
7850cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
786e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
787e78576d6SMark F. Adams     if (jj > 0) {
7880cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7890cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7900cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7912e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7922e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7932e68589bSMark F. Adams       PetscInt       *fids;
79465d7b583SSatish Balay       PetscReal      *data;
7950cbbd2e1SMark F. Adams       /* count agg */
7960cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7970cbbd2e1SMark F. Adams 
7980cbbd2e1SMark F. Adams       /* get block */
7992e68589bSMark F. Adams       ierr = PetscMalloc((Mdata*N)*sizeof(PetscScalar), &qqc);CHKERRQ(ierr);
8002e68589bSMark F. Adams       ierr = PetscMalloc((M*N)*sizeof(PetscScalar), &qqr);CHKERRQ(ierr);
8012e68589bSMark F. Adams       ierr = PetscMalloc(N*sizeof(PetscScalar), &TAU);CHKERRQ(ierr);
8022e68589bSMark F. Adams       ierr = PetscMalloc(LWORK*sizeof(PetscScalar), &WORK);CHKERRQ(ierr);
8032e68589bSMark F. Adams       ierr = PetscMalloc(M*sizeof(PetscInt), &fids);CHKERRQ(ierr);
8042e68589bSMark F. Adams 
8052e68589bSMark F. Adams       aggID = 0;
806e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
807e78576d6SMark F. Adams       while (pos) {
808e78576d6SMark F. Adams         PetscInt gid1;
809ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
810e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
811e78576d6SMark F. Adams 
8120cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
8130cbbd2e1SMark F. Adams         else {
8140cbbd2e1SMark F. Adams           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
81571959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
8160cbbd2e1SMark F. Adams         }
8172e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
81865d7b583SSatish Balay         data = &data_in[flid*bs];
8192e68589bSMark F. Adams         for (kk = ii = 0; ii < bs; ii++) {
8202e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
8210cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8220cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8232e68589bSMark F. Adams           }
8242e68589bSMark F. Adams         }
825519f805aSKarl Rupp #if defined(OUT_AGGS)
826b2a4f308SMark F. Adams         if (llev==1) {
8279057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8280cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
829c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
830f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8310cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
832b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
833b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8349057884aSMark F. Adams         }
8359057884aSMark F. Adams #endif
8362e68589bSMark F. Adams         /* set fine IDs */
8372e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8382e68589bSMark F. Adams 
8392e68589bSMark F. Adams         aggID++;
8400cbbd2e1SMark F. Adams       }
8412e68589bSMark F. Adams 
8422e68589bSMark F. Adams       /* pad with zeros */
8432e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8442e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8452e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8462e68589bSMark F. Adams         }
8472e68589bSMark F. Adams       }
8482e68589bSMark F. Adams 
8492e68589bSMark F. Adams       ndone += aggID;
8502e68589bSMark F. Adams       /* QR */
85184df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
852f75e95b9SBarry Smith       PetscStackCall("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
85384df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
854d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8552e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8562e68589bSMark F. Adams       {
8572e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8582e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8592e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
86071959b99SBarry Smith            if (data[jj*out_data_stride + ii] != 1.e300) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != 1.e300");
8610cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8620cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8632e68589bSMark F. Adams           }
8642e68589bSMark F. Adams         }
8652e68589bSMark F. Adams       }
8662e68589bSMark F. Adams 
8672e68589bSMark F. Adams       /* get Q - row oriented */
868f75e95b9SBarry Smith       PetscStackCall("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
869c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8702e68589bSMark F. Adams 
8712e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8722e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8732e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8742e68589bSMark F. Adams         }
8752e68589bSMark F. Adams       }
8762e68589bSMark F. Adams 
8772e68589bSMark F. Adams       /* add diagonal block of P0 */
878c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
879c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
880c8b0795cSMark F. Adams       }
8812e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8822e68589bSMark F. Adams 
8832e68589bSMark F. Adams       ierr = PetscFree(qqc);CHKERRQ(ierr);
8842e68589bSMark F. Adams       ierr = PetscFree(qqr);CHKERRQ(ierr);
8852e68589bSMark F. Adams       ierr = PetscFree(TAU);CHKERRQ(ierr);
8862e68589bSMark F. Adams       ierr = PetscFree(WORK);CHKERRQ(ierr);
8872e68589bSMark F. Adams       ierr = PetscFree(fids);CHKERRQ(ierr);
888b43b03e9SMark F. Adams       clid++;
8890cbbd2e1SMark F. Adams     } /* coarse agg */
8900cbbd2e1SMark F. Adams   } /* for all fine nodes */
8910cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8920cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8932e68589bSMark F. Adams 
8943b4367a7SBarry Smith /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8952e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
8963b4367a7SBarry Smith /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8973b4367a7SBarry 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); */
8982e68589bSMark F. Adams 
899519f805aSKarl Rupp #if defined(OUT_AGGS)
900b2a4f308SMark F. Adams   if (llev==1) fclose(file);
9019057884aSMark F. Adams #endif
9020cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
9032e68589bSMark F. Adams   PetscFunctionReturn(0);
9042e68589bSMark F. Adams }
9052e68589bSMark F. Adams 
9062e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
9072e68589bSMark F. Adams /*
908c8b0795cSMark F. Adams    PCGAMGgraph_AGG
9092e68589bSMark F. Adams 
9102e68589bSMark F. Adams   Input Parameter:
9112e68589bSMark F. Adams    . pc - this
9122e68589bSMark F. Adams    . Amat - matrix on this fine level
9132e68589bSMark F. Adams   Output Parameter:
914c8b0795cSMark F. Adams    . a_Gmat -
9152e68589bSMark F. Adams */
9162e68589bSMark F. Adams #undef __FUNCT__
917c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG"
9182fa5cd67SKarl Rupp PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat)
919c8b0795cSMark F. Adams {
920c8b0795cSMark F. Adams   PetscErrorCode            ierr;
921c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
922c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
923c8b0795cSMark F. Adams   const PetscInt            verbose      = pc_gamg->verbose;
924c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
925c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
926c5df96a5SBarry Smith   PetscMPIInt               rank,size;
927e0940f08SMark F. Adams   Mat                       Gmat;
9283b4367a7SBarry Smith   MPI_Comm                  comm;
92967744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
930c8b0795cSMark F. Adams 
931c8b0795cSMark F. Adams   PetscFunctionBegin;
9323b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
9330cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9340cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9350cbbd2e1SMark F. Adams #endif
9363b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9373b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
938c8b0795cSMark F. Adams 
93967744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
940c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9410cbbd2e1SMark F. Adams 
9422d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
9432d7fac45SMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
944c8b0795cSMark F. Adams 
945e0940f08SMark F. Adams   *a_Gmat = Gmat;
946c8b0795cSMark F. Adams 
9470cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9480cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9490cbbd2e1SMark F. Adams #endif
950c8b0795cSMark F. Adams   PetscFunctionReturn(0);
951c8b0795cSMark F. Adams }
952c8b0795cSMark F. Adams 
953c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
954c8b0795cSMark F. Adams /*
955b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
956c8b0795cSMark F. Adams 
957c8b0795cSMark F. Adams   Input Parameter:
958e0940f08SMark F. Adams    . a_pc - this
959e0940f08SMark F. Adams   Input/Output Parameter:
9600cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
961c8b0795cSMark F. Adams   Output Parameter:
9620cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
963c8b0795cSMark F. Adams */
964c8b0795cSMark F. Adams #undef __FUNCT__
965b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
9662fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
967c8b0795cSMark F. Adams {
968c8b0795cSMark F. Adams   PetscErrorCode ierr;
969e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
970c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
971c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9720cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9739f3f12c8SMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
9740cbbd2e1SMark F. Adams   IS             perm;
975c8b0795cSMark F. Adams   PetscInt       Ii,nloc,bs,n,m;
976c8b0795cSMark F. Adams   PetscInt       *permute;
977c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
978b43b03e9SMark F. Adams   MatCoarsen     crs;
9793b4367a7SBarry Smith   MPI_Comm       comm;
980c5df96a5SBarry Smith   PetscMPIInt    rank,size;
981c8b0795cSMark F. Adams 
982c8b0795cSMark F. Adams   PetscFunctionBegin;
9830cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9840cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9850cbbd2e1SMark F. Adams #endif
9863b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9873b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9883b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
989e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
99071959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
99171959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
992c8b0795cSMark F. Adams   nloc = n/bs;
993c8b0795cSMark F. Adams 
994e0940f08SMark F. Adams   if (pc_gamg_agg->square_graph) {
9953b4367a7SBarry Smith     if (verbose > 1) PetscPrintf(comm,"[%d]%s square graph\n",rank,__FUNCT__);
996878e152fSMark F. Adams     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
997806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
9989f3f12c8SMark F. Adams     if (verbose > 2) {
9993b4367a7SBarry Smith       ierr = PetscPrintf(comm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr);
1000f10ff945SMark F. Adams       /* check for symetry */
1001f10ff945SMark F. Adams       if (verbose > 4) {
1002f10ff945SMark F. Adams 
1003f10ff945SMark F. Adams       }
10049f3f12c8SMark F. Adams     }
1005806fa848SBarry Smith   } else Gmat2 = Gmat1;
1006c8b0795cSMark F. Adams 
1007c8b0795cSMark F. Adams   /* get MIS aggs */
1008c8b0795cSMark F. Adams   /* randomize */
1009c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
1010c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
1011c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1012c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
1013c8b0795cSMark F. Adams     permute[Ii]   = Ii;
1014c8b0795cSMark F. Adams   }
1015c8b0795cSMark F. Adams   srand(1); /* make deterministic */
1016c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1017c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
1018c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1019c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1020c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1021c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1022c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1023c8b0795cSMark F. Adams     }
1024c8b0795cSMark F. Adams   }
1025c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
1026c8b0795cSMark F. Adams 
10273b4367a7SBarry Smith   if (verbose > 1) PetscPrintf(comm,"[%d]%s coarsen graph\n",rank,__FUNCT__);
10289f3f12c8SMark F. Adams 
1029806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
10300cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10310cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1032b43b03e9SMark F. Adams #endif
10333b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
10349057884aSMark F. Adams   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
10359057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1036b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1037b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1038b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
1039b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1040b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
10410cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1042b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1043b43b03e9SMark F. Adams 
1044c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1045c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10460cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10470cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1048b43b03e9SMark F. Adams #endif
10493b4367a7SBarry Smith   if (verbose > 2) PetscPrintf(comm,"[%d]%s coarsen graph done\n",rank,__FUNCT__);
10509f3f12c8SMark F. Adams 
1051c8b0795cSMark F. Adams   /* smooth aggs */
1052e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10530cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10540cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1055c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1056e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
105741b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10583b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1059806fa848SBarry Smith   } else {
10600cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1061f10ff945SMark F. Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
106241b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10630cbbd2e1SMark F. Adams     if (mat) {
10640cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10650cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10660cbbd2e1SMark F. Adams     }
10670cbbd2e1SMark F. Adams   }
10680cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
10690cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
10700cbbd2e1SMark F. Adams #endif
10713b4367a7SBarry Smith   if (verbose > 2) PetscPrintf(comm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__);
1072c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1073c8b0795cSMark F. Adams }
1074c8b0795cSMark F. Adams 
1075c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1076c8b0795cSMark F. Adams /*
10770cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1078c8b0795cSMark F. Adams 
1079c8b0795cSMark F. Adams  Input Parameter:
1080c8b0795cSMark F. Adams  . pc - this
1081c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1082c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10830cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1084c8b0795cSMark F. Adams  Output Parameter:
1085c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1086c8b0795cSMark F. Adams  */
1087c8b0795cSMark F. Adams #undef __FUNCT__
10880cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
10892fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10902e68589bSMark F. Adams {
10912e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10922e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
10932e68589bSMark F. Adams   const PetscInt verbose   = pc_gamg->verbose;
1094c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10952e68589bSMark F. Adams   PetscErrorCode ierr;
1096c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1097c8b0795cSMark F. Adams   Mat            Prol;
1098c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10993b4367a7SBarry Smith   MPI_Comm       comm;
11000cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1101c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1102c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
11032e68589bSMark F. Adams 
11042e68589bSMark F. Adams   PetscFunctionBegin;
11053b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
11060cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11070cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11080cbbd2e1SMark F. Adams #endif
11093b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
11103b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
11112e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1112c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
111371959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
111471959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
11152e68589bSMark F. Adams 
11162e68589bSMark F. Adams   /* get 'nLocalSelected' */
11170cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1118e78576d6SMark F. Adams     PetscBool ise;
11190cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1120e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1121e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
11222e68589bSMark F. Adams   }
11232e68589bSMark F. Adams 
11242e68589bSMark F. Adams   /* create prolongator, create P matrix */
11253b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1126806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1127a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1128a2f3521dSMark F. Adams   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
11290298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
11300298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1131a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1132a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
11330298fd71SBarry Smith   /* data_cols, NULL, data_cols, NULL, */
1134a2f3521dSMark F. Adams   /* &Prol); */
11352e68589bSMark F. Adams 
11362e68589bSMark F. Adams   /* can get all points "removed" */
1137c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1138c8b0795cSMark F. Adams   if (ii==0) {
11393b4367a7SBarry Smith     if (verbose > 0) PetscPrintf(comm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__);
11402e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
11410298fd71SBarry Smith     *a_P_out = NULL;  /* out */
11422e68589bSMark F. Adams     PetscFunctionReturn(0);
11432e68589bSMark F. Adams   }
11443b4367a7SBarry Smith   if (verbose > 0) PetscPrintf(comm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs);
1145c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
11460cbbd2e1SMark F. Adams 
114771959b99SBarry 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);
1148c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
114971959b99SBarry 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);
11502e68589bSMark F. Adams 
11512e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11530cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11542e68589bSMark F. Adams #endif
1155c5df96a5SBarry Smith   if (size > 1) { /*  */
11562e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
11572e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr);
11582e68589bSMark F. Adams     for (jj = 0; jj < data_cols; jj++) {
1159c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1160a2f3521dSMark F. Adams         PetscInt        ii,stride;
1161c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11622fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11632fa5cd67SKarl Rupp 
1164806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1165a2f3521dSMark F. Adams 
11662e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1167a2f3521dSMark F. Adams           ierr    = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr);
1168a2f3521dSMark F. Adams           nbnodes = bs*stride;
11692e68589bSMark F. Adams         }
1170a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11712fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11722e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11732e68589bSMark F. Adams       }
11742e68589bSMark F. Adams     }
11752e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1176806fa848SBarry Smith   } else {
1177c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1178c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11792e68589bSMark F. Adams   }
11802e68589bSMark F. Adams 
11812e68589bSMark F. Adams   /* get P0 */
1182c5df96a5SBarry Smith   if (size > 1) {
11832e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1184a2f3521dSMark F. Adams     PetscInt  stride;
11852e68589bSMark F. Adams 
11862e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr);
11872e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1188806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1189a2f3521dSMark F. Adams     ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1190a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11912e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1192a2f3521dSMark F. Adams 
119371959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11942e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1195806fa848SBarry Smith   } else {
11962e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
11972e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11982e68589bSMark F. Adams   }
11990cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12000cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
12010cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
12022e68589bSMark F. Adams #endif
1203c8b0795cSMark F. Adams   {
12040298fd71SBarry Smith     PetscReal *data_out = NULL;
1205806fa848SBarry Smith     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1206c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12072fa5cd67SKarl Rupp 
1208c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
1209c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1210c8b0795cSMark F. Adams     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1211c8b0795cSMark F. Adams   }
12120cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12130cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1214c8b0795cSMark F. Adams #endif
1215c5df96a5SBarry Smith   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
12162e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
12172e68589bSMark F. Adams 
1218c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
12190cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12200cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
12210cbbd2e1SMark F. Adams #endif
1222c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1223c8b0795cSMark F. Adams }
1224c8b0795cSMark F. Adams 
1225c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1226c8b0795cSMark F. Adams /*
12270cbbd2e1SMark F. Adams    PCGAMGOptprol_AGG
1228c8b0795cSMark F. Adams 
1229c8b0795cSMark F. Adams   Input Parameter:
1230c8b0795cSMark F. Adams    . pc - this
1231c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1232c8b0795cSMark F. Adams  In/Output Parameter:
1233c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1234c8b0795cSMark F. Adams */
1235c8b0795cSMark F. Adams #undef __FUNCT__
12360cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG"
12372fa5cd67SKarl Rupp PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P)
1238c8b0795cSMark F. Adams {
1239c8b0795cSMark F. Adams   PetscErrorCode ierr;
1240c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1241c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1242c8b0795cSMark F. Adams   const PetscInt verbose      = pc_gamg->verbose;
1243c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1244c8b0795cSMark F. Adams   PetscInt       jj;
1245c5df96a5SBarry Smith   PetscMPIInt    rank,size;
1246c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
12473b4367a7SBarry Smith   MPI_Comm       comm;
1248c8b0795cSMark F. Adams 
1249c8b0795cSMark F. Adams   PetscFunctionBegin;
12503b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
12510cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12520cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
12530cbbd2e1SMark F. Adams #endif
12543b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12553b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1256c8b0795cSMark F. Adams 
12572e68589bSMark F. Adams   /* smooth P0 */
1258c8b0795cSMark F. Adams   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12592e68589bSMark F. Adams     Mat       tMat;
12602e68589bSMark F. Adams     Vec       diag;
12612e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12630cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12642e68589bSMark F. Adams #endif
12652e68589bSMark F. Adams     if (jj == 0) {
12662e68589bSMark F. Adams       KSP eksp;
12672e68589bSMark F. Adams       Vec bb, xx;
1268da509fc8SJed Brown       PC  epc;
12692e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr);
12702e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr);
12712e68589bSMark F. Adams       {
12722e68589bSMark F. Adams         PetscRandom rctx;
12733b4367a7SBarry Smith         ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
12742e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12752e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12762e68589bSMark F. Adams         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
12772e68589bSMark F. Adams       }
1278e696ed0bSMark F. Adams 
1279e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1280e696ed0bSMark F. Adams       {
1281e696ed0bSMark F. Adams         PetscInt    Istart,Iend,ncols,jj,Ii;
1282e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1283e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1284e696ed0bSMark F. Adams         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1285e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1286e696ed0bSMark F. Adams           if (ncols <= 1) {
1287e696ed0bSMark F. Adams             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1288e696ed0bSMark F. Adams           }
1289e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1290e696ed0bSMark F. Adams         }
1291e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1292e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1293e696ed0bSMark F. Adams       }
1294e696ed0bSMark F. Adams 
12953b4367a7SBarry Smith       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1296806fa848SBarry Smith       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12972e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1298da509fc8SJed Brown       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1299c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1300c2436cacSMark F. Adams       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1301c2436cacSMark F. Adams 
1302c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1303806fa848SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13042e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
13052e68589bSMark F. Adams 
1306da509fc8SJed Brown       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1307da509fc8SJed Brown       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1308c2436cacSMark F. Adams 
13092e68589bSMark F. Adams       /* solve - keep stuff out of logging */
13102e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
13112e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
13122e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
13132e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
13142e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
13152e68589bSMark F. Adams 
13162e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
13179f3f12c8SMark F. Adams       if (verbose > 0) {
13183b4367a7SBarry Smith         PetscPrintf(comm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
13192e68589bSMark F. Adams                     __FUNCT__,emax,emin,PCJACOBI);
13202e68589bSMark F. Adams       }
13212e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
13222e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
13232e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
13242e68589bSMark F. Adams 
13252e68589bSMark F. Adams       if (pc_gamg->emax_id == -1) {
13263bf036e2SBarry Smith         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
132771959b99SBarry Smith         if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1");
13282e68589bSMark F. Adams       }
13292e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
13302e68589bSMark F. Adams     }
13312e68589bSMark F. Adams 
13322e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
13332e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
13342e68589bSMark F. Adams     ierr  = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr);
13352e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
13362e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
13372e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
13382e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1339b4ec6429SMark F. Adams     alpha = -1.4/emax;
13402e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
13412e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
13422e68589bSMark F. Adams     Prol  = tMat;
13430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13440cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
13452e68589bSMark F. Adams #endif
13462e68589bSMark F. Adams   }
13470cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13480cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
13490cbbd2e1SMark F. Adams #endif
1350c8b0795cSMark F. Adams   *a_P = Prol;
1351c8b0795cSMark F. Adams   PetscFunctionReturn(0);
13522e68589bSMark F. Adams }
13532e68589bSMark F. Adams 
1354c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1355c8b0795cSMark F. Adams /*
1356a2f3521dSMark F. Adams    PCGAMGKKTProl_AGG
1357a2f3521dSMark F. Adams 
1358a2f3521dSMark F. Adams   Input Parameter:
1359a2f3521dSMark F. Adams    . pc - this
1360a2f3521dSMark F. Adams    . Prol11 - matrix on this fine level
1361a2f3521dSMark F. Adams    . A21 - matrix on this fine level
1362a2f3521dSMark F. Adams  In/Output Parameter:
1363a2f3521dSMark F. Adams    . a_P22 - prolongation operator to the next level
1364a2f3521dSMark F. Adams */
1365a2f3521dSMark F. Adams #undef __FUNCT__
1366a2f3521dSMark F. Adams #define __FUNCT__ "PCGAMGKKTProl_AGG"
13672fa5cd67SKarl Rupp PetscErrorCode PCGAMGKKTProl_AGG(PC pc,const Mat Prol11,const Mat A21,Mat *a_P22)
1368a2f3521dSMark F. Adams {
1369a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1370a2f3521dSMark F. Adams   PC_MG            *mg      = (PC_MG*)pc->data;
1371a2f3521dSMark F. Adams   PC_GAMG          *pc_gamg = (PC_GAMG*)mg->innerctx;
1372a2f3521dSMark F. Adams   const PetscInt   verbose  = pc_gamg->verbose;
1373a2f3521dSMark F. Adams   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1374c5df96a5SBarry Smith   PetscMPIInt      rank,size;
1375a2f3521dSMark F. Adams   Mat              Prol22,Tmat,Gmat;
13763b4367a7SBarry Smith   MPI_Comm         comm;
1377a2f3521dSMark F. Adams   PetscCoarsenData *agg_lists;
1378a2f3521dSMark F. Adams 
1379a2f3521dSMark F. Adams   PetscFunctionBegin;
13803b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1381a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1382a2f3521dSMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1383a2f3521dSMark F. Adams #endif
13843b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
13853b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1386a2f3521dSMark F. Adams 
1387a2f3521dSMark F. Adams   /* form C graph */
1388a2f3521dSMark F. Adams   ierr = MatMatMult(A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);CHKERRQ(ierr);
1389a2f3521dSMark F. Adams   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat);CHKERRQ(ierr);
1390a2f3521dSMark F. Adams   ierr = MatDestroy(&Tmat);CHKERRQ(ierr);
1391a2f3521dSMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose);CHKERRQ(ierr);
1392a2f3521dSMark F. Adams 
1393a2f3521dSMark F. Adams   /* coarsen constraints */
1394a2f3521dSMark F. Adams   {
1395a2f3521dSMark F. Adams     MatCoarsen crs;
13963b4367a7SBarry Smith     ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
1397a2f3521dSMark F. Adams     ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
1398a2f3521dSMark F. Adams     ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
1399a2f3521dSMark F. Adams     ierr = MatCoarsenSetVerbose(crs, verbose);CHKERRQ(ierr);
1400a2f3521dSMark F. Adams     ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1401a2f3521dSMark F. Adams     ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1402a2f3521dSMark F. Adams     ierr = MatCoarsenGetData(crs, &agg_lists);CHKERRQ(ierr);
1403a2f3521dSMark F. Adams     ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1404a2f3521dSMark F. Adams   }
1405a2f3521dSMark F. Adams 
1406a2f3521dSMark F. Adams   /* form simple prolongation 'Prol22' */
1407a2f3521dSMark F. Adams   {
1408a2f3521dSMark F. Adams     PetscInt    ii,mm,clid,my0,nloc,nLocalSelected;
1409a2f3521dSMark F. Adams     PetscScalar val = 1.0;
1410a2f3521dSMark F. Adams     /* get 'nLocalSelected' */
1411a2f3521dSMark F. Adams     ierr = MatGetLocalSize(Gmat, &nloc, &ii);CHKERRQ(ierr);
1412a2f3521dSMark F. Adams     for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1413a2f3521dSMark F. Adams       PetscBool ise;
1414a2f3521dSMark F. Adams       /* filter out singletons 0 or 1? */
1415a2f3521dSMark F. Adams       ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1416a2f3521dSMark F. Adams       if (!ise) nLocalSelected++;
1417a2f3521dSMark F. Adams     }
1418a2f3521dSMark F. Adams 
14193b4367a7SBarry Smith     ierr = MatCreate(comm,&Prol22);CHKERRQ(ierr);
1420806fa848SBarry Smith     ierr = MatSetSizes(Prol22,nloc, nLocalSelected,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1421a2f3521dSMark F. Adams     ierr = MatSetType(Prol22, MATAIJ);CHKERRQ(ierr);
14220298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(Prol22,1,NULL);CHKERRQ(ierr);
14230298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(Prol22,1,NULL,1,NULL);CHKERRQ(ierr);
14243b4367a7SBarry Smith     /* ierr = MatCreateAIJ(comm, */
1425a2f3521dSMark F. Adams     /*                      nloc, nLocalSelected, */
1426a2f3521dSMark F. Adams     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
14270298fd71SBarry Smith     /*                      1, NULL, 1, NULL, */
1428a2f3521dSMark F. Adams     /*                      &Prol22); */
1429a2f3521dSMark F. Adams 
1430a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(Prol22, &my0, &ii);CHKERRQ(ierr);
1431a2f3521dSMark F. Adams     nloc = ii - my0;
1432a2f3521dSMark F. Adams 
1433a2f3521dSMark F. Adams     /* make aggregates */
1434a2f3521dSMark F. Adams     for (mm = clid = 0; mm < nloc; mm++) {
1435a2f3521dSMark F. Adams       ierr = PetscCDSizeAt(agg_lists, mm, &ii);CHKERRQ(ierr);
1436a2f3521dSMark F. Adams       if (ii > 0) {
1437a2f3521dSMark F. Adams         PetscInt   asz=ii,cgid=my0+clid,rids[1000];
1438a2f3521dSMark F. Adams         PetscCDPos pos;
1439c666adbfSMark F. Adams         if (asz>1000) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1440a2f3521dSMark F. Adams         ii   = 0;
1441a2f3521dSMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1442a2f3521dSMark F. Adams         while (pos) {
1443a2f3521dSMark F. Adams           PetscInt gid1;
1444a2f3521dSMark F. Adams           ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
1445a2f3521dSMark F. Adams           ierr = PetscCDGetNextPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1446a2f3521dSMark F. Adams 
1447a2f3521dSMark F. Adams           rids[ii++] = gid1;
1448a2f3521dSMark F. Adams         }
144971959b99SBarry Smith         if (ii != asz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D != asz %D",ii,asz);
1450a2f3521dSMark F. Adams         /* add diagonal block of P0 */
1451a2f3521dSMark F. Adams         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES);CHKERRQ(ierr);
1452a2f3521dSMark F. Adams 
1453a2f3521dSMark F. Adams         clid++;
1454a2f3521dSMark F. Adams       } /* coarse agg */
1455a2f3521dSMark F. Adams     } /* for all fine nodes */
1456a2f3521dSMark F. Adams     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1457a2f3521dSMark F. Adams     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1458a2f3521dSMark F. Adams   }
1459a2f3521dSMark F. Adams 
1460a2f3521dSMark F. Adams   /* clean up */
1461a2f3521dSMark F. Adams   ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
1462a2f3521dSMark F. Adams   ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
1463a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1464a2f3521dSMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1465a2f3521dSMark F. Adams #endif
1466a2f3521dSMark F. Adams   *a_P22 = Prol22;
1467a2f3521dSMark F. Adams   PetscFunctionReturn(0);
1468a2f3521dSMark F. Adams }
1469a2f3521dSMark F. Adams 
1470a2f3521dSMark F. Adams /* -------------------------------------------------------------------------- */
1471a2f3521dSMark F. Adams /*
1472c8b0795cSMark F. Adams    PCCreateGAMG_AGG
14732e68589bSMark F. Adams 
1474c8b0795cSMark F. Adams   Input Parameter:
1475c8b0795cSMark F. Adams    . pc -
1476c8b0795cSMark F. Adams */
1477c8b0795cSMark F. Adams #undef __FUNCT__
1478c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1479c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1480c8b0795cSMark F. Adams {
1481c8b0795cSMark F. Adams   PetscErrorCode ierr;
1482c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1483c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1484c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
14852e68589bSMark F. Adams 
1486c8b0795cSMark F. Adams   PetscFunctionBegin;
1487ef820a72SMark F. Adams   if (pc_gamg->subctx) {
1488ef820a72SMark F. Adams     /* call base class */
1489ef820a72SMark F. Adams     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
1490ef820a72SMark F. Adams   }
1491ef820a72SMark F. Adams 
1492c8b0795cSMark F. Adams   /* create sub context for SA */
1493c8b0795cSMark F. Adams   ierr            = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr);
1494c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1495c8b0795cSMark F. Adams 
1496*1ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1497c8b0795cSMark F. Adams   pc->ops->destroy        = PCDestroy_AGG;
1498c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1499c8b0795cSMark F. Adams 
1500c8b0795cSMark F. Adams   /* set internal function pointers */
1501*1ab5ffc9SJed Brown   pc_gamg->ops->graph       = PCGAMGgraph_AGG;
1502*1ab5ffc9SJed Brown   pc_gamg->ops->coarsen     = PCGAMGCoarsen_AGG;
1503*1ab5ffc9SJed Brown   pc_gamg->ops->prolongator = PCGAMGProlongator_AGG;
1504*1ab5ffc9SJed Brown   pc_gamg->ops->optprol     = PCGAMGOptprol_AGG;
1505*1ab5ffc9SJed Brown   pc_gamg->ops->formkktprol = PCGAMGKKTProl_AGG;
1506c8b0795cSMark F. Adams 
1507*1ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
1508c8b0795cSMark F. Adams 
1509bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
15102e68589bSMark F. Adams   PetscFunctionReturn(0);
15112e68589bSMark F. Adams }
1512