xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 8b83055f287c8deb1a16325ae3faee24d0648b7a)
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__
2059b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG"
2069b8ffb57SJed Brown PetscErrorCode PCDestroy_GAMG_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;
2112e68589bSMark F. Adams 
2122e68589bSMark F. Adams   PetscFunctionBegin;
2139b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
2142e68589bSMark F. Adams   PetscFunctionReturn(0);
2152e68589bSMark F. Adams }
2162e68589bSMark F. Adams 
2172e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2182e68589bSMark F. Adams /*
2192e68589bSMark F. Adams    PCSetCoordinates_AGG
220302f38e8SMark F. Adams      - collective
2212e68589bSMark F. Adams 
2222e68589bSMark F. Adams    Input Parameter:
2232e68589bSMark F. Adams    . pc - the preconditioner context
224a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
225302f38e8SMark F. Adams    . a_nloc - number of vertices local
226302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2272e68589bSMark F. Adams */
2281e6b0712SBarry Smith 
2292e68589bSMark F. Adams #undef __FUNCT__
2302e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2311e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2322e68589bSMark F. Adams {
2332e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2342e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2352e68589bSMark F. Adams   PetscErrorCode ierr;
23669344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
237a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2382e68589bSMark F. Adams 
2392e68589bSMark F. Adams   PetscFunctionBegin;
240a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
241a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
242302f38e8SMark F. Adams   nloc = a_nloc;
2432e68589bSMark F. Adams 
2442e68589bSMark F. Adams   /* SA: null space vectors */
24569344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
24669344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
247a2f3521dSMark F. Adams   else if (coords) {
248c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
24969344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
25069344418SMark F. Adams     if (ndm != ndf) {
251c666adbfSMark 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);
252a2f3521dSMark F. Adams     }
253806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
25471959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
25571959b99SBarry 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);
256c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2572e68589bSMark F. Adams 
2582e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2592e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2602e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
261302f38e8SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
2622e68589bSMark F. Adams   }
2632e68589bSMark F. Adams   /* copy data in - column oriented */
2642e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
26569344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
26669344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
267c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2682e68589bSMark F. Adams     else {
26969344418SMark F. Adams       /* translational modes */
2702fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2712fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
27269344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2732e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2742fa5cd67SKarl Rupp         }
2752fa5cd67SKarl Rupp       }
27669344418SMark F. Adams 
27769344418SMark F. Adams       /* rotational modes */
2782e68589bSMark F. Adams       if (coords) {
27969344418SMark F. Adams         if (ndm == 2) {
2802e68589bSMark F. Adams           data   += 2*M;
2812e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2822e68589bSMark F. Adams           data[1] =  coords[2*kk];
283806fa848SBarry Smith         } else {
2842e68589bSMark F. Adams           data   += 3*M;
2852e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2862e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2872e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2882e68589bSMark F. Adams         }
2892e68589bSMark F. Adams       }
2902e68589bSMark F. Adams     }
2912e68589bSMark F. Adams   }
2922e68589bSMark F. Adams 
2932e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2942e68589bSMark F. Adams   PetscFunctionReturn(0);
2952e68589bSMark F. Adams }
2962e68589bSMark F. Adams 
297b43b03e9SMark F. Adams typedef PetscInt NState;
298b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
299b43b03e9SMark F. Adams static const NState DELETED =-1;
300b43b03e9SMark F. Adams static const NState REMOVED =-3;
301b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
302b43b03e9SMark F. Adams 
303c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
304c8b0795cSMark F. Adams /*
305b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
306b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
307c8b0795cSMark F. Adams 
308c8b0795cSMark F. Adams    Input Parameter:
309c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
310c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
311c8b0795cSMark F. Adams    Input/Output Parameter:
3120cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
313c8b0795cSMark F. Adams */
314c8b0795cSMark F. Adams #undef __FUNCT__
315c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3160cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
3170cbbd2e1SMark F. Adams                                  const Mat Gmat_1,  /* base graph */
3180cbbd2e1SMark F. Adams                                  /* const IS selected_2, [nselected local] selected vertices */
3191147fc2aSKarl Rupp                                  PetscCoarsenData *aggs_2)  /* [nselected local] global ID of aggregate */
320c8b0795cSMark F. Adams {
321c8b0795cSMark F. Adams   PetscErrorCode ierr;
322c8b0795cSMark F. Adams   PetscBool      isMPI;
323c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
3243b4367a7SBarry Smith   MPI_Comm       comm;
325c5df96a5SBarry Smith   PetscMPIInt    rank,size;
3260cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
327c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
328c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3290cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3300cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
331c8b0795cSMark F. Adams   NState         *lid_state;
3320cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
333c8b0795cSMark F. Adams 
334c8b0795cSMark F. Adams   PetscFunctionBegin;
3353b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
3363b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
3373b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
338c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
339c8b0795cSMark F. Adams 
3400cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
341c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
342c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3433b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
344c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
345c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
346c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
347c8b0795cSMark F. Adams   }
348c8b0795cSMark F. Adams 
349c8b0795cSMark F. Adams   /* get submatrices */
350251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
351c8b0795cSMark F. Adams   if (isMPI) {
352c8b0795cSMark F. Adams     /* grab matrix objects */
353c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
354c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
355c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
356c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
357c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
358c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
359c8b0795cSMark F. Adams 
360c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
361c8b0795cSMark F. Adams     matB_1->compressedrow.check = PETSC_TRUE;
3622fa5cd67SKarl Rupp 
363806fa848SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
364c8b0795cSMark F. Adams 
365c8b0795cSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &lid_cprowID_1);CHKERRQ(ierr);
3660cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
367c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
368c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
369c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
370c8b0795cSMark F. Adams     }
371806fa848SBarry Smith   } else {
372c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
373c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3740298fd71SBarry Smith     lid_cprowID_1 = NULL;
375c8b0795cSMark F. Adams   }
37671959b99SBarry Smith   if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
37771959b99SBarry Smith   if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
37871959b99SBarry Smith   if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
37971959b99SBarry Smith   if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
380c8b0795cSMark F. Adams 
381c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
382c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(NState), &lid_state);CHKERRQ(ierr);
3830cbbd2e1SMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscScalar), &lid_parent_gid);CHKERRQ(ierr);
384c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3850cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
386c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
387c8b0795cSMark F. Adams   }
3880cbbd2e1SMark F. Adams 
3890cbbd2e1SMark F. Adams   /* set lid_state */
3900cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
39141b27cdeSMark F. Adams     PetscCDPos pos;
392e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
393e78576d6SMark F. Adams     if (pos) {
394e78576d6SMark F. Adams       PetscInt gid1;
3952fa5cd67SKarl Rupp 
39671959b99SBarry Smith       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
39771959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3980cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
399b43b03e9SMark F. Adams     }
400b43b03e9SMark F. Adams   }
4010cbbd2e1SMark F. Adams 
4020cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
403c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
404c8b0795cSMark F. Adams     NState state = lid_state[lid];
405c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
40641b27cdeSMark F. Adams       PetscCDPos pos;
407e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
408e78576d6SMark F. Adams       while (pos) {
409e78576d6SMark F. Adams         PetscInt gid1;
410ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
411e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
412e78576d6SMark F. Adams 
4132fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
414c8b0795cSMark F. Adams       }
4150cbbd2e1SMark F. Adams     }
4160cbbd2e1SMark F. Adams   }
4170cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
418c8b0795cSMark F. Adams   if (isMPI) {
419c8b0795cSMark F. Adams     Vec tempVec;
420c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
421c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
422c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
423c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
424c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
425c8b0795cSMark F. Adams     }
426c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
427c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
428806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
429806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
430c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
431c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
432806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
433806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
434c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4350cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4360cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4370cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4380cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4390cbbd2e1SMark F. Adams     }
4400cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4410cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4420cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
443806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
444806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4450cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4460cbbd2e1SMark F. Adams 
447c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
448c8b0795cSMark F. Adams   } /* ismpi */
449c8b0795cSMark F. Adams 
450c8b0795cSMark F. Adams   /* doit */
451c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
452c8b0795cSMark F. Adams     NState state = lid_state[lid];
4530cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4540cbbd2e1SMark F. Adams       /* steal locals */
455c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
456c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
457c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4580cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
459c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4600cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4610cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4620cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4630cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
4640298fd71SBarry Smith             PetscCDPos pos,last=NULL;
465c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
466e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
467e78576d6SMark F. Adams             while (pos) {
468e78576d6SMark F. Adams               PetscInt gid;
469ffc955d6SMark F. Adams               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4700cbbd2e1SMark F. Adams               if (gid == gidj) {
47171959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
47241b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
47341b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4740cbbd2e1SMark F. Adams                 hav  = 1;
475c8b0795cSMark F. Adams                 break;
476806fa848SBarry Smith               } else last = pos;
477e78576d6SMark F. Adams 
478e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
479c8b0795cSMark F. Adams             }
480c8b0795cSMark F. Adams             if (hav!=1) {
48171959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
482c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
483c8b0795cSMark F. Adams             }
484806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4852fa5cd67SKarl 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.");
48641b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
487c8b0795cSMark F. Adams           }
488c8b0795cSMark F. Adams         }
4890cbbd2e1SMark F. Adams       } /* local neighbors */
490806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4910cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
492c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
493c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
494c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
495c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
496c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
497c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
498c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
499c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
5000cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
5010cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
5020cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
5030298fd71SBarry Smith               PetscCDPos pos,last=NULL;
5040cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
505e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
506e78576d6SMark F. Adams               while (pos) {
507e78576d6SMark F. Adams                 PetscInt gid;
508ffc955d6SMark F. Adams                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
5090cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
5100cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
51171959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
51241b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
5130cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
5140cbbd2e1SMark F. Adams                   hav = 1;
515c8b0795cSMark F. Adams                   break;
516806fa848SBarry Smith                 } else last = pos;
517e78576d6SMark F. Adams 
518e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
519c8b0795cSMark F. Adams               }
520c8b0795cSMark F. Adams               if (hav!=1) {
521c666adbfSMark F. Adams                 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
522c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
523c8b0795cSMark F. Adams               }
524806fa848SBarry Smith             } else {
5250cbbd2e1SMark F. Adams               /* ghosts remove this later */
5260cbbd2e1SMark F. Adams             }
527c8b0795cSMark F. Adams           }
528c8b0795cSMark F. Adams         }
529c8b0795cSMark F. Adams       }
530c8b0795cSMark F. Adams     } /* selected/deleted */
531c8b0795cSMark F. Adams   } /* node loop */
532c8b0795cSMark F. Adams 
533c8b0795cSMark F. Adams   if (isMPI) {
5340cbbd2e1SMark F. Adams     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
5350cbbd2e1SMark F. Adams     Vec           tempVec,ghostgids2,ghostparents2;
5360cbbd2e1SMark F. Adams     PetscInt      cpid,nghost_2;
5370cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
538c8b0795cSMark F. Adams 
5390cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
540c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5410cbbd2e1SMark F. Adams 
5420cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
543c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5440cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
545c8b0795cSMark F. Adams     }
546c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
547c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5480cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
549806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
550806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5510cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5520cbbd2e1SMark F. Adams 
5530cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5540cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5550cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5560cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5570cbbd2e1SMark F. Adams     }
5580cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5590cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5600cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
561806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
562806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5630cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5640cbbd2e1SMark F. Adams 
565c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
566c8b0795cSMark F. Adams 
5670cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
568f10ff945SMark F. Adams     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5690cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5700cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5710cbbd2e1SMark F. Adams       if (state==DELETED) {
5720cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5730cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5740cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5750cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5760cbbd2e1SMark F. Adams           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5770cbbd2e1SMark F. Adams         }
5780cbbd2e1SMark F. Adams       }
5790cbbd2e1SMark F. Adams     }
580c8b0795cSMark F. Adams 
5810cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
582c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
583c8b0795cSMark F. Adams       NState state = lid_state[lid];
584c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
5850298fd71SBarry Smith         PetscCDPos pos,last=NULL;
586c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
587e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
588e78576d6SMark F. Adams         while (pos) {
589e78576d6SMark F. Adams           PetscInt gid;
590ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
591e78576d6SMark F. Adams 
5920cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5930cbbd2e1SMark F. Adams             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5940cbbd2e1SMark F. Adams             if (cpid != -1) {
5950cbbd2e1SMark F. Adams               /* a moved ghost - */
5960cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
59741b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
598806fa848SBarry Smith             } else last = pos;
599806fa848SBarry Smith           } else last = pos;
600e78576d6SMark F. Adams 
601e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
602c8b0795cSMark F. Adams         } /* loop over list of deleted */
603c8b0795cSMark F. Adams       } /* selected */
604c8b0795cSMark F. Adams     }
6050cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
606c8b0795cSMark F. Adams 
6070cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
6080cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
6090cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6100cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
6110cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6120cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
61341b27cdeSMark F. Adams         PetscCDPos pos;
6140cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
615e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
616e78576d6SMark F. Adams         while (pos) {
617e78576d6SMark F. Adams           PetscInt gidj;
618ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
619e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
620e78576d6SMark F. Adams 
6210cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
622c8b0795cSMark F. Adams         }
623c8b0795cSMark F. Adams         if (hav != 1) {
624ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
62541b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
626c8b0795cSMark F. Adams         }
627c8b0795cSMark F. Adams       }
628c8b0795cSMark F. Adams     }
629c8b0795cSMark F. Adams 
6300cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6310cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6320cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6330cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
634c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6350cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6360cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6370cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
638c8b0795cSMark F. Adams   }
639c8b0795cSMark F. Adams 
6400cbbd2e1SMark F. Adams   ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr);
641c8b0795cSMark F. Adams   ierr = PetscFree(lid_state);CHKERRQ(ierr);
642c8b0795cSMark F. Adams   PetscFunctionReturn(0);
643c8b0795cSMark F. Adams }
6442e68589bSMark F. Adams 
6452e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6462e68589bSMark F. Adams /*
647a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
648a2f3521dSMark F. Adams       Looks in Mat for near null space.
649a2f3521dSMark F. Adams       Does not work for Stokes
6502e68589bSMark F. Adams 
6512e68589bSMark F. Adams   Input Parameter:
6522e68589bSMark F. Adams    . pc -
653a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6542e68589bSMark F. Adams */
6552e68589bSMark F. Adams #undef __FUNCT__
6562e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
657b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6582e68589bSMark F. Adams {
6592e68589bSMark F. Adams   PetscErrorCode ierr;
660b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
661b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
662b8cd405aSMark F. Adams   MatNullSpace   mnull;
663b8cd405aSMark F. Adams 
6642e68589bSMark F. Adams   PetscFunctionBegin;
665b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
666b8cd405aSMark F. Adams   if (!mnull) {
667a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6689d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
66971959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
67071959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6710298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
672806fa848SBarry Smith   } else {
673b8cd405aSMark F. Adams     PetscReal *nullvec;
674b8cd405aSMark F. Adams     PetscBool has_const;
675b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
676b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6770298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
678b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
6798caf3d72SBarry Smith     ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
680b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
681b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
682b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
683b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
684b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
685b8cd405aSMark F. Adams     }
686b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
687b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6882fa5cd67SKarl Rupp 
689a2f3521dSMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); /* this does not work for Stokes */
6902fa5cd67SKarl Rupp 
691b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
692b8cd405aSMark F. Adams   }
6932e68589bSMark F. Adams   PetscFunctionReturn(0);
6942e68589bSMark F. Adams }
6952e68589bSMark F. Adams 
6962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6972e68589bSMark F. Adams /*
6982e68589bSMark F. Adams  formProl0
6992e68589bSMark F. Adams 
7002e68589bSMark F. Adams    Input Parameter:
7010cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
7022e68589bSMark F. Adams    . bs - block size
7030cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7040cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7052e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
7062e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7072e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7082e68589bSMark F. Adams   Output Parameter:
7092e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7102e68589bSMark F. Adams    . a_Prol - prolongation operator
7112e68589bSMark F. Adams */
7122e68589bSMark F. Adams #undef __FUNCT__
7132e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7140cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */
7150cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
7160cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
7170cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
7180cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7190cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7200cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7212e68589bSMark F. Adams                                 PetscReal **a_data_out,
7221147fc2aSKarl Rupp                                 Mat a_Prol) /* prolongation operator (output)*/
7232e68589bSMark F. Adams {
7242e68589bSMark F. Adams   PetscErrorCode ierr;
7250cbbd2e1SMark F. Adams   PetscInt       Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7263b4367a7SBarry Smith   MPI_Comm       comm;
727c5df96a5SBarry Smith   PetscMPIInt    rank, size;
7282e68589bSMark F. Adams   PetscReal      *out_data;
72941b27cdeSMark F. Adams   PetscCDPos     pos;
7300cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7310cbbd2e1SMark F. Adams 
732797e13b7SMark F. Adams /* #define OUT_AGGS */
733519f805aSKarl Rupp #if defined(OUT_AGGS)
7340298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7359057884aSMark F. Adams #endif
7362e68589bSMark F. Adams 
7372e68589bSMark F. Adams   PetscFunctionBegin;
7383b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7393b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7403b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
7412e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
74271959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
74371959b99SBarry 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);
7440cbbd2e1SMark F. Adams   Iend   /= bs;
7450cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7462e68589bSMark F. Adams 
747f10ff945SMark F. Adams   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7480cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7490cbbd2e1SMark F. Adams     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7502e68589bSMark F. Adams   }
7512e68589bSMark F. Adams 
752519f805aSKarl Rupp #if defined(OUT_AGGS)
753c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7542fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7550cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7560cbbd2e1SMark F. Adams #endif
7570cbbd2e1SMark F. Adams 
7580cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7590cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
760e78576d6SMark F. Adams     PetscBool ise;
761e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
762e78576d6SMark F. Adams     if (!ise) nSelected++;
7630cbbd2e1SMark F. Adams   }
7640cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
76571959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
76671959b99SBarry 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);
7670cbbd2e1SMark F. Adams 
7682e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7690cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7702fa5cd67SKarl Rupp 
7710cbbd2e1SMark F. Adams   ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr);
7722fa5cd67SKarl Rupp   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=1.e300;
7730cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7742e68589bSMark F. Adams 
7752e68589bSMark F. Adams   /* find points and set prolongation */
776c8b0795cSMark F. Adams   minsz = 100;
7772e68589bSMark F. Adams   ndone = 0;
7780cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
779e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
780e78576d6SMark F. Adams     if (jj > 0) {
7810cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7820cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7830cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7842e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7852e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7862e68589bSMark F. Adams       PetscInt       *fids;
78765d7b583SSatish Balay       PetscReal      *data;
7880cbbd2e1SMark F. Adams       /* count agg */
7890cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7900cbbd2e1SMark F. Adams 
7910cbbd2e1SMark F. Adams       /* get block */
7922e68589bSMark F. Adams       ierr = PetscMalloc((Mdata*N)*sizeof(PetscScalar), &qqc);CHKERRQ(ierr);
7932e68589bSMark F. Adams       ierr = PetscMalloc((M*N)*sizeof(PetscScalar), &qqr);CHKERRQ(ierr);
7942e68589bSMark F. Adams       ierr = PetscMalloc(N*sizeof(PetscScalar), &TAU);CHKERRQ(ierr);
7952e68589bSMark F. Adams       ierr = PetscMalloc(LWORK*sizeof(PetscScalar), &WORK);CHKERRQ(ierr);
7962e68589bSMark F. Adams       ierr = PetscMalloc(M*sizeof(PetscInt), &fids);CHKERRQ(ierr);
7972e68589bSMark F. Adams 
7982e68589bSMark F. Adams       aggID = 0;
799e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
800e78576d6SMark F. Adams       while (pos) {
801e78576d6SMark F. Adams         PetscInt gid1;
802ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
803e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
804e78576d6SMark F. Adams 
8050cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
8060cbbd2e1SMark F. Adams         else {
8070cbbd2e1SMark F. Adams           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
80871959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
8090cbbd2e1SMark F. Adams         }
8102e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
81165d7b583SSatish Balay         data = &data_in[flid*bs];
8122e68589bSMark F. Adams         for (kk = ii = 0; ii < bs; ii++) {
8132e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
8140cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8150cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8162e68589bSMark F. Adams           }
8172e68589bSMark F. Adams         }
818519f805aSKarl Rupp #if defined(OUT_AGGS)
819b2a4f308SMark F. Adams         if (llev==1) {
8209057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8210cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
822c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
823f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8240cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
825b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
826b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8279057884aSMark F. Adams         }
8289057884aSMark F. Adams #endif
8292e68589bSMark F. Adams         /* set fine IDs */
8302e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8312e68589bSMark F. Adams 
8322e68589bSMark F. Adams         aggID++;
8330cbbd2e1SMark F. Adams       }
8342e68589bSMark F. Adams 
8352e68589bSMark F. Adams       /* pad with zeros */
8362e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8372e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8382e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8392e68589bSMark F. Adams         }
8402e68589bSMark F. Adams       }
8412e68589bSMark F. Adams 
8422e68589bSMark F. Adams       ndone += aggID;
8432e68589bSMark F. Adams       /* QR */
84484df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
845*8b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
84684df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
847d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8482e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8492e68589bSMark F. Adams       {
8502e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8512e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8522e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
85371959b99SBarry 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");
8540cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8550cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8562e68589bSMark F. Adams           }
8572e68589bSMark F. Adams         }
8582e68589bSMark F. Adams       }
8592e68589bSMark F. Adams 
8602e68589bSMark F. Adams       /* get Q - row oriented */
861*8b83055fSJed Brown       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
862c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8632e68589bSMark F. Adams 
8642e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8652e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8662e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8672e68589bSMark F. Adams         }
8682e68589bSMark F. Adams       }
8692e68589bSMark F. Adams 
8702e68589bSMark F. Adams       /* add diagonal block of P0 */
871c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
872c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
873c8b0795cSMark F. Adams       }
8742e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8752e68589bSMark F. Adams 
8762e68589bSMark F. Adams       ierr = PetscFree(qqc);CHKERRQ(ierr);
8772e68589bSMark F. Adams       ierr = PetscFree(qqr);CHKERRQ(ierr);
8782e68589bSMark F. Adams       ierr = PetscFree(TAU);CHKERRQ(ierr);
8792e68589bSMark F. Adams       ierr = PetscFree(WORK);CHKERRQ(ierr);
8802e68589bSMark F. Adams       ierr = PetscFree(fids);CHKERRQ(ierr);
881b43b03e9SMark F. Adams       clid++;
8820cbbd2e1SMark F. Adams     } /* coarse agg */
8830cbbd2e1SMark F. Adams   } /* for all fine nodes */
8840cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8850cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8862e68589bSMark F. Adams 
8873b4367a7SBarry Smith /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8882e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
8893b4367a7SBarry Smith /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8903b4367a7SBarry 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); */
8912e68589bSMark F. Adams 
892519f805aSKarl Rupp #if defined(OUT_AGGS)
893b2a4f308SMark F. Adams   if (llev==1) fclose(file);
8949057884aSMark F. Adams #endif
8950cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
8962e68589bSMark F. Adams   PetscFunctionReturn(0);
8972e68589bSMark F. Adams }
8982e68589bSMark F. Adams 
8992e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
9002e68589bSMark F. Adams /*
901c8b0795cSMark F. Adams    PCGAMGgraph_AGG
9022e68589bSMark F. Adams 
9032e68589bSMark F. Adams   Input Parameter:
9042e68589bSMark F. Adams    . pc - this
9052e68589bSMark F. Adams    . Amat - matrix on this fine level
9062e68589bSMark F. Adams   Output Parameter:
907c8b0795cSMark F. Adams    . a_Gmat -
9082e68589bSMark F. Adams */
9092e68589bSMark F. Adams #undef __FUNCT__
910c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG"
9112fa5cd67SKarl Rupp PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat)
912c8b0795cSMark F. Adams {
913c8b0795cSMark F. Adams   PetscErrorCode            ierr;
914c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
915c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
916c8b0795cSMark F. Adams   const PetscInt            verbose      = pc_gamg->verbose;
917c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
918c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
919c5df96a5SBarry Smith   PetscMPIInt               rank,size;
920e0940f08SMark F. Adams   Mat                       Gmat;
9213b4367a7SBarry Smith   MPI_Comm                  comm;
92267744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
923c8b0795cSMark F. Adams 
924c8b0795cSMark F. Adams   PetscFunctionBegin;
9253b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
9260cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9270cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9280cbbd2e1SMark F. Adams #endif
9293b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9303b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
931c8b0795cSMark F. Adams 
93267744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
933c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9340cbbd2e1SMark F. Adams 
9352d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
9362d7fac45SMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
937c8b0795cSMark F. Adams 
938e0940f08SMark F. Adams   *a_Gmat = Gmat;
939c8b0795cSMark F. Adams 
9400cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9410cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9420cbbd2e1SMark F. Adams #endif
943c8b0795cSMark F. Adams   PetscFunctionReturn(0);
944c8b0795cSMark F. Adams }
945c8b0795cSMark F. Adams 
946c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
947c8b0795cSMark F. Adams /*
948b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
949c8b0795cSMark F. Adams 
950c8b0795cSMark F. Adams   Input Parameter:
951e0940f08SMark F. Adams    . a_pc - this
952e0940f08SMark F. Adams   Input/Output Parameter:
9530cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
954c8b0795cSMark F. Adams   Output Parameter:
9550cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
956c8b0795cSMark F. Adams */
957c8b0795cSMark F. Adams #undef __FUNCT__
958b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
9592fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
960c8b0795cSMark F. Adams {
961c8b0795cSMark F. Adams   PetscErrorCode ierr;
962e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
963c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
964c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9650cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9669f3f12c8SMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
9670cbbd2e1SMark F. Adams   IS             perm;
968c8b0795cSMark F. Adams   PetscInt       Ii,nloc,bs,n,m;
969c8b0795cSMark F. Adams   PetscInt       *permute;
970c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
971b43b03e9SMark F. Adams   MatCoarsen     crs;
9723b4367a7SBarry Smith   MPI_Comm       comm;
973c5df96a5SBarry Smith   PetscMPIInt    rank,size;
974c8b0795cSMark F. Adams 
975c8b0795cSMark F. Adams   PetscFunctionBegin;
9760cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9770cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9780cbbd2e1SMark F. Adams #endif
9793b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9803b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9813b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
982e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
98371959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
98471959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
985c8b0795cSMark F. Adams   nloc = n/bs;
986c8b0795cSMark F. Adams 
987e0940f08SMark F. Adams   if (pc_gamg_agg->square_graph) {
9883b4367a7SBarry Smith     if (verbose > 1) PetscPrintf(comm,"[%d]%s square graph\n",rank,__FUNCT__);
989878e152fSMark F. Adams     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
990806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
9919f3f12c8SMark F. Adams     if (verbose > 2) {
9923b4367a7SBarry Smith       ierr = PetscPrintf(comm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr);
993f10ff945SMark F. Adams       /* check for symetry */
994f10ff945SMark F. Adams       if (verbose > 4) {
995f10ff945SMark F. Adams 
996f10ff945SMark F. Adams       }
9979f3f12c8SMark F. Adams     }
998806fa848SBarry Smith   } else Gmat2 = Gmat1;
999c8b0795cSMark F. Adams 
1000c8b0795cSMark F. Adams   /* get MIS aggs */
1001c8b0795cSMark F. Adams   /* randomize */
1002c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
1003c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
1004c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1005c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
1006c8b0795cSMark F. Adams     permute[Ii]   = Ii;
1007c8b0795cSMark F. Adams   }
1008c8b0795cSMark F. Adams   srand(1); /* make deterministic */
1009c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1010c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
1011c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1012c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1013c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1014c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1015c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1016c8b0795cSMark F. Adams     }
1017c8b0795cSMark F. Adams   }
1018c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
1019c8b0795cSMark F. Adams 
10203b4367a7SBarry Smith   if (verbose > 1) PetscPrintf(comm,"[%d]%s coarsen graph\n",rank,__FUNCT__);
10219f3f12c8SMark F. Adams 
1022806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
10230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10240cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1025b43b03e9SMark F. Adams #endif
10263b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
10279057884aSMark F. Adams   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
10289057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1029b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1030b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1031b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
1032b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1033b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
10340cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1035b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1036b43b03e9SMark F. Adams 
1037c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1038c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10400cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1041b43b03e9SMark F. Adams #endif
10423b4367a7SBarry Smith   if (verbose > 2) PetscPrintf(comm,"[%d]%s coarsen graph done\n",rank,__FUNCT__);
10439f3f12c8SMark F. Adams 
1044c8b0795cSMark F. Adams   /* smooth aggs */
1045e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10460cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10470cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1048c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1049e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
105041b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10513b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1052806fa848SBarry Smith   } else {
10530cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1054f10ff945SMark F. Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
105541b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10560cbbd2e1SMark F. Adams     if (mat) {
10570cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10580cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10590cbbd2e1SMark F. Adams     }
10600cbbd2e1SMark F. Adams   }
10610cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
10620cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
10630cbbd2e1SMark F. Adams #endif
10643b4367a7SBarry Smith   if (verbose > 2) PetscPrintf(comm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__);
1065c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1066c8b0795cSMark F. Adams }
1067c8b0795cSMark F. Adams 
1068c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1069c8b0795cSMark F. Adams /*
10700cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1071c8b0795cSMark F. Adams 
1072c8b0795cSMark F. Adams  Input Parameter:
1073c8b0795cSMark F. Adams  . pc - this
1074c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1075c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10760cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1077c8b0795cSMark F. Adams  Output Parameter:
1078c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1079c8b0795cSMark F. Adams  */
1080c8b0795cSMark F. Adams #undef __FUNCT__
10810cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
10822fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10832e68589bSMark F. Adams {
10842e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10852e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
10862e68589bSMark F. Adams   const PetscInt verbose   = pc_gamg->verbose;
1087c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10882e68589bSMark F. Adams   PetscErrorCode ierr;
1089c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1090c8b0795cSMark F. Adams   Mat            Prol;
1091c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10923b4367a7SBarry Smith   MPI_Comm       comm;
10930cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1094c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1095c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
10962e68589bSMark F. Adams 
10972e68589bSMark F. Adams   PetscFunctionBegin;
10983b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
10990cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11000cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11010cbbd2e1SMark F. Adams #endif
11023b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
11033b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
11042e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1105c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
110671959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
110771959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
11082e68589bSMark F. Adams 
11092e68589bSMark F. Adams   /* get 'nLocalSelected' */
11100cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1111e78576d6SMark F. Adams     PetscBool ise;
11120cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1113e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1114e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
11152e68589bSMark F. Adams   }
11162e68589bSMark F. Adams 
11172e68589bSMark F. Adams   /* create prolongator, create P matrix */
11183b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1119806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1120a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1121a2f3521dSMark F. Adams   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
11220298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
11230298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1124a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1125a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
11260298fd71SBarry Smith   /* data_cols, NULL, data_cols, NULL, */
1127a2f3521dSMark F. Adams   /* &Prol); */
11282e68589bSMark F. Adams 
11292e68589bSMark F. Adams   /* can get all points "removed" */
1130c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1131c8b0795cSMark F. Adams   if (ii==0) {
11323b4367a7SBarry Smith     if (verbose > 0) PetscPrintf(comm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__);
11332e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
11340298fd71SBarry Smith     *a_P_out = NULL;  /* out */
11352e68589bSMark F. Adams     PetscFunctionReturn(0);
11362e68589bSMark F. Adams   }
11373b4367a7SBarry Smith   if (verbose > 0) PetscPrintf(comm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs);
1138c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
11390cbbd2e1SMark F. Adams 
114071959b99SBarry 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);
1141c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
114271959b99SBarry 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);
11432e68589bSMark F. Adams 
11442e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11460cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11472e68589bSMark F. Adams #endif
1148c5df96a5SBarry Smith   if (size > 1) { /*  */
11492e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
11502e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr);
11512e68589bSMark F. Adams     for (jj = 0; jj < data_cols; jj++) {
1152c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1153a2f3521dSMark F. Adams         PetscInt        ii,stride;
1154c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11552fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11562fa5cd67SKarl Rupp 
1157806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1158a2f3521dSMark F. Adams 
11592e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1160a2f3521dSMark F. Adams           ierr    = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr);
1161a2f3521dSMark F. Adams           nbnodes = bs*stride;
11622e68589bSMark F. Adams         }
1163a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11642fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11652e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11662e68589bSMark F. Adams       }
11672e68589bSMark F. Adams     }
11682e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1169806fa848SBarry Smith   } else {
1170c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1171c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11722e68589bSMark F. Adams   }
11732e68589bSMark F. Adams 
11742e68589bSMark F. Adams   /* get P0 */
1175c5df96a5SBarry Smith   if (size > 1) {
11762e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1177a2f3521dSMark F. Adams     PetscInt  stride;
11782e68589bSMark F. Adams 
11792e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr);
11802e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1181806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1182a2f3521dSMark F. Adams     ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1183a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11842e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1185a2f3521dSMark F. Adams 
118671959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11872e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1188806fa848SBarry Smith   } else {
11892e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
11902e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11912e68589bSMark F. Adams   }
11920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11930cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11940cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
11952e68589bSMark F. Adams #endif
1196c8b0795cSMark F. Adams   {
11970298fd71SBarry Smith     PetscReal *data_out = NULL;
1198806fa848SBarry Smith     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1199c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12002fa5cd67SKarl Rupp 
1201c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
1202c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1203c8b0795cSMark F. Adams     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1204c8b0795cSMark F. Adams   }
12050cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12060cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1207c8b0795cSMark F. Adams #endif
1208c5df96a5SBarry Smith   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
12092e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
12102e68589bSMark F. Adams 
1211c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
12120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12130cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
12140cbbd2e1SMark F. Adams #endif
1215c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1216c8b0795cSMark F. Adams }
1217c8b0795cSMark F. Adams 
1218c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1219c8b0795cSMark F. Adams /*
12200cbbd2e1SMark F. Adams    PCGAMGOptprol_AGG
1221c8b0795cSMark F. Adams 
1222c8b0795cSMark F. Adams   Input Parameter:
1223c8b0795cSMark F. Adams    . pc - this
1224c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1225c8b0795cSMark F. Adams  In/Output Parameter:
1226c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1227c8b0795cSMark F. Adams */
1228c8b0795cSMark F. Adams #undef __FUNCT__
12290cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG"
12302fa5cd67SKarl Rupp PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P)
1231c8b0795cSMark F. Adams {
1232c8b0795cSMark F. Adams   PetscErrorCode ierr;
1233c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1234c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1235c8b0795cSMark F. Adams   const PetscInt verbose      = pc_gamg->verbose;
1236c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1237c8b0795cSMark F. Adams   PetscInt       jj;
1238c5df96a5SBarry Smith   PetscMPIInt    rank,size;
1239c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
12403b4367a7SBarry Smith   MPI_Comm       comm;
1241c8b0795cSMark F. Adams 
1242c8b0795cSMark F. Adams   PetscFunctionBegin;
12433b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
12440cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12450cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
12460cbbd2e1SMark F. Adams #endif
12473b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12483b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1249c8b0795cSMark F. Adams 
12502e68589bSMark F. Adams   /* smooth P0 */
1251c8b0795cSMark F. Adams   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12522e68589bSMark F. Adams     Mat       tMat;
12532e68589bSMark F. Adams     Vec       diag;
12542e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12550cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12560cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12572e68589bSMark F. Adams #endif
12582e68589bSMark F. Adams     if (jj == 0) {
12592e68589bSMark F. Adams       KSP eksp;
12602e68589bSMark F. Adams       Vec bb, xx;
1261da509fc8SJed Brown       PC  epc;
12622e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr);
12632e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr);
12642e68589bSMark F. Adams       {
12652e68589bSMark F. Adams         PetscRandom rctx;
12663b4367a7SBarry Smith         ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
12672e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12682e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12692e68589bSMark F. Adams         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
12702e68589bSMark F. Adams       }
1271e696ed0bSMark F. Adams 
1272e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1273e696ed0bSMark F. Adams       {
1274e696ed0bSMark F. Adams         PetscInt    Istart,Iend,ncols,jj,Ii;
1275e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1276e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1277e696ed0bSMark F. Adams         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1278e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1279e696ed0bSMark F. Adams           if (ncols <= 1) {
1280e696ed0bSMark F. Adams             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1281e696ed0bSMark F. Adams           }
1282e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1283e696ed0bSMark F. Adams         }
1284e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1285e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1286e696ed0bSMark F. Adams       }
1287e696ed0bSMark F. Adams 
12883b4367a7SBarry Smith       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1289806fa848SBarry Smith       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12902e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1291da509fc8SJed Brown       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1292c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1293c2436cacSMark F. Adams       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1294c2436cacSMark F. Adams 
1295c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1296806fa848SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
12972e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
12982e68589bSMark F. Adams 
1299da509fc8SJed Brown       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1300da509fc8SJed Brown       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1301c2436cacSMark F. Adams 
13022e68589bSMark F. Adams       /* solve - keep stuff out of logging */
13032e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
13042e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
13052e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
13062e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
13072e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
13082e68589bSMark F. Adams 
13092e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
13109f3f12c8SMark F. Adams       if (verbose > 0) {
13113b4367a7SBarry Smith         PetscPrintf(comm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
13122e68589bSMark F. Adams                     __FUNCT__,emax,emin,PCJACOBI);
13132e68589bSMark F. Adams       }
13142e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
13152e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
13162e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
13172e68589bSMark F. Adams 
13182e68589bSMark F. Adams       if (pc_gamg->emax_id == -1) {
13193bf036e2SBarry Smith         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
132071959b99SBarry Smith         if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1");
13212e68589bSMark F. Adams       }
13222e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
13232e68589bSMark F. Adams     }
13242e68589bSMark F. Adams 
13252e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
13262e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
13272e68589bSMark F. Adams     ierr  = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr);
13282e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
13292e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
13302e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
13312e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1332b4ec6429SMark F. Adams     alpha = -1.4/emax;
13332e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
13342e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
13352e68589bSMark F. Adams     Prol  = tMat;
13360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13370cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
13382e68589bSMark F. Adams #endif
13392e68589bSMark F. Adams   }
13400cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13410cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
13420cbbd2e1SMark F. Adams #endif
1343c8b0795cSMark F. Adams   *a_P = Prol;
1344c8b0795cSMark F. Adams   PetscFunctionReturn(0);
13452e68589bSMark F. Adams }
13462e68589bSMark F. Adams 
1347c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1348c8b0795cSMark F. Adams /*
1349a2f3521dSMark F. Adams    PCGAMGKKTProl_AGG
1350a2f3521dSMark F. Adams 
1351a2f3521dSMark F. Adams   Input Parameter:
1352a2f3521dSMark F. Adams    . pc - this
1353a2f3521dSMark F. Adams    . Prol11 - matrix on this fine level
1354a2f3521dSMark F. Adams    . A21 - matrix on this fine level
1355a2f3521dSMark F. Adams  In/Output Parameter:
1356a2f3521dSMark F. Adams    . a_P22 - prolongation operator to the next level
1357a2f3521dSMark F. Adams */
1358a2f3521dSMark F. Adams #undef __FUNCT__
1359a2f3521dSMark F. Adams #define __FUNCT__ "PCGAMGKKTProl_AGG"
13602fa5cd67SKarl Rupp PetscErrorCode PCGAMGKKTProl_AGG(PC pc,const Mat Prol11,const Mat A21,Mat *a_P22)
1361a2f3521dSMark F. Adams {
1362a2f3521dSMark F. Adams   PetscErrorCode   ierr;
1363a2f3521dSMark F. Adams   PC_MG            *mg      = (PC_MG*)pc->data;
1364a2f3521dSMark F. Adams   PC_GAMG          *pc_gamg = (PC_GAMG*)mg->innerctx;
1365a2f3521dSMark F. Adams   const PetscInt   verbose  = pc_gamg->verbose;
1366a2f3521dSMark F. Adams   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1367c5df96a5SBarry Smith   PetscMPIInt      rank,size;
1368a2f3521dSMark F. Adams   Mat              Prol22,Tmat,Gmat;
13693b4367a7SBarry Smith   MPI_Comm         comm;
1370a2f3521dSMark F. Adams   PetscCoarsenData *agg_lists;
1371a2f3521dSMark F. Adams 
1372a2f3521dSMark F. Adams   PetscFunctionBegin;
13733b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1374a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1375a2f3521dSMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1376a2f3521dSMark F. Adams #endif
13773b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
13783b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1379a2f3521dSMark F. Adams 
1380a2f3521dSMark F. Adams   /* form C graph */
1381a2f3521dSMark F. Adams   ierr = MatMatMult(A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);CHKERRQ(ierr);
1382a2f3521dSMark F. Adams   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat);CHKERRQ(ierr);
1383a2f3521dSMark F. Adams   ierr = MatDestroy(&Tmat);CHKERRQ(ierr);
1384a2f3521dSMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose);CHKERRQ(ierr);
1385a2f3521dSMark F. Adams 
1386a2f3521dSMark F. Adams   /* coarsen constraints */
1387a2f3521dSMark F. Adams   {
1388a2f3521dSMark F. Adams     MatCoarsen crs;
13893b4367a7SBarry Smith     ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
1390a2f3521dSMark F. Adams     ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
1391a2f3521dSMark F. Adams     ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
1392a2f3521dSMark F. Adams     ierr = MatCoarsenSetVerbose(crs, verbose);CHKERRQ(ierr);
1393a2f3521dSMark F. Adams     ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1394a2f3521dSMark F. Adams     ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1395a2f3521dSMark F. Adams     ierr = MatCoarsenGetData(crs, &agg_lists);CHKERRQ(ierr);
1396a2f3521dSMark F. Adams     ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1397a2f3521dSMark F. Adams   }
1398a2f3521dSMark F. Adams 
1399a2f3521dSMark F. Adams   /* form simple prolongation 'Prol22' */
1400a2f3521dSMark F. Adams   {
1401a2f3521dSMark F. Adams     PetscInt    ii,mm,clid,my0,nloc,nLocalSelected;
1402a2f3521dSMark F. Adams     PetscScalar val = 1.0;
1403a2f3521dSMark F. Adams     /* get 'nLocalSelected' */
1404a2f3521dSMark F. Adams     ierr = MatGetLocalSize(Gmat, &nloc, &ii);CHKERRQ(ierr);
1405a2f3521dSMark F. Adams     for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1406a2f3521dSMark F. Adams       PetscBool ise;
1407a2f3521dSMark F. Adams       /* filter out singletons 0 or 1? */
1408a2f3521dSMark F. Adams       ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1409a2f3521dSMark F. Adams       if (!ise) nLocalSelected++;
1410a2f3521dSMark F. Adams     }
1411a2f3521dSMark F. Adams 
14123b4367a7SBarry Smith     ierr = MatCreate(comm,&Prol22);CHKERRQ(ierr);
1413806fa848SBarry Smith     ierr = MatSetSizes(Prol22,nloc, nLocalSelected,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1414a2f3521dSMark F. Adams     ierr = MatSetType(Prol22, MATAIJ);CHKERRQ(ierr);
14150298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(Prol22,1,NULL);CHKERRQ(ierr);
14160298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(Prol22,1,NULL,1,NULL);CHKERRQ(ierr);
14173b4367a7SBarry Smith     /* ierr = MatCreateAIJ(comm, */
1418a2f3521dSMark F. Adams     /*                      nloc, nLocalSelected, */
1419a2f3521dSMark F. Adams     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
14200298fd71SBarry Smith     /*                      1, NULL, 1, NULL, */
1421a2f3521dSMark F. Adams     /*                      &Prol22); */
1422a2f3521dSMark F. Adams 
1423a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(Prol22, &my0, &ii);CHKERRQ(ierr);
1424a2f3521dSMark F. Adams     nloc = ii - my0;
1425a2f3521dSMark F. Adams 
1426a2f3521dSMark F. Adams     /* make aggregates */
1427a2f3521dSMark F. Adams     for (mm = clid = 0; mm < nloc; mm++) {
1428a2f3521dSMark F. Adams       ierr = PetscCDSizeAt(agg_lists, mm, &ii);CHKERRQ(ierr);
1429a2f3521dSMark F. Adams       if (ii > 0) {
1430a2f3521dSMark F. Adams         PetscInt   asz=ii,cgid=my0+clid,rids[1000];
1431a2f3521dSMark F. Adams         PetscCDPos pos;
1432c666adbfSMark F. Adams         if (asz>1000) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1433a2f3521dSMark F. Adams         ii   = 0;
1434a2f3521dSMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1435a2f3521dSMark F. Adams         while (pos) {
1436a2f3521dSMark F. Adams           PetscInt gid1;
1437a2f3521dSMark F. Adams           ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
1438a2f3521dSMark F. Adams           ierr = PetscCDGetNextPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1439a2f3521dSMark F. Adams 
1440a2f3521dSMark F. Adams           rids[ii++] = gid1;
1441a2f3521dSMark F. Adams         }
144271959b99SBarry Smith         if (ii != asz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D != asz %D",ii,asz);
1443a2f3521dSMark F. Adams         /* add diagonal block of P0 */
1444a2f3521dSMark F. Adams         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES);CHKERRQ(ierr);
1445a2f3521dSMark F. Adams 
1446a2f3521dSMark F. Adams         clid++;
1447a2f3521dSMark F. Adams       } /* coarse agg */
1448a2f3521dSMark F. Adams     } /* for all fine nodes */
1449a2f3521dSMark F. Adams     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1450a2f3521dSMark F. Adams     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1451a2f3521dSMark F. Adams   }
1452a2f3521dSMark F. Adams 
1453a2f3521dSMark F. Adams   /* clean up */
1454a2f3521dSMark F. Adams   ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
1455a2f3521dSMark F. Adams   ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
1456a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1457a2f3521dSMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1458a2f3521dSMark F. Adams #endif
1459a2f3521dSMark F. Adams   *a_P22 = Prol22;
1460a2f3521dSMark F. Adams   PetscFunctionReturn(0);
1461a2f3521dSMark F. Adams }
1462a2f3521dSMark F. Adams 
1463a2f3521dSMark F. Adams /* -------------------------------------------------------------------------- */
1464a2f3521dSMark F. Adams /*
1465c8b0795cSMark F. Adams    PCCreateGAMG_AGG
14662e68589bSMark F. Adams 
1467c8b0795cSMark F. Adams   Input Parameter:
1468c8b0795cSMark F. Adams    . pc -
1469c8b0795cSMark F. Adams */
1470c8b0795cSMark F. Adams #undef __FUNCT__
1471c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1472c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1473c8b0795cSMark F. Adams {
1474c8b0795cSMark F. Adams   PetscErrorCode ierr;
1475c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1476c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1477c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
14782e68589bSMark F. Adams 
1479c8b0795cSMark F. Adams   PetscFunctionBegin;
1480c8b0795cSMark F. Adams   /* create sub context for SA */
1481c8b0795cSMark F. Adams   ierr            = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr);
1482c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1483c8b0795cSMark F. Adams 
14841ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
14859b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1486c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1487c8b0795cSMark F. Adams 
1488c8b0795cSMark F. Adams   /* set internal function pointers */
14891ab5ffc9SJed Brown   pc_gamg->ops->graph       = PCGAMGgraph_AGG;
14901ab5ffc9SJed Brown   pc_gamg->ops->coarsen     = PCGAMGCoarsen_AGG;
14911ab5ffc9SJed Brown   pc_gamg->ops->prolongator = PCGAMGProlongator_AGG;
14921ab5ffc9SJed Brown   pc_gamg->ops->optprol     = PCGAMGOptprol_AGG;
14931ab5ffc9SJed Brown   pc_gamg->ops->formkktprol = PCGAMGKKTProl_AGG;
1494c8b0795cSMark F. Adams 
14951ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
1496c8b0795cSMark F. Adams 
1497bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
14982e68589bSMark F. Adams   PetscFunctionReturn(0);
14992e68589bSMark F. Adams }
1500