xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
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 <assert.h>
92e68589bSMark F. Adams #include <petscblaslapack.h>
102e68589bSMark F. Adams 
112e68589bSMark F. Adams typedef struct {
12c8b0795cSMark F. Adams   PetscInt  nsmooths;
13c8b0795cSMark F. Adams   PetscBool sym_graph;
14ef4ad70eSMark F. Adams   PetscBool square_graph;
152e68589bSMark F. Adams } PC_GAMG_AGG;
162e68589bSMark F. Adams 
172e68589bSMark F. Adams #undef __FUNCT__
182e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths"
192e68589bSMark F. Adams /*@
202e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
212e68589bSMark F. Adams 
222e68589bSMark F. Adams    Not Collective on PC
232e68589bSMark F. Adams 
242e68589bSMark F. Adams    Input Parameters:
252e68589bSMark F. Adams .  pc - the preconditioner context
262e68589bSMark F. Adams 
272e68589bSMark F. Adams    Options Database Key:
282e68589bSMark F. Adams .  -pc_gamg_agg_nsmooths
292e68589bSMark F. Adams 
302e68589bSMark F. Adams    Level: intermediate
312e68589bSMark F. Adams 
322e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
332e68589bSMark F. Adams 
342e68589bSMark F. Adams .seealso: ()
352e68589bSMark F. Adams @*/
362e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
372e68589bSMark F. Adams {
382e68589bSMark F. Adams   PetscErrorCode ierr;
392e68589bSMark F. Adams 
402e68589bSMark F. Adams   PetscFunctionBegin;
412e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
422e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
432e68589bSMark F. Adams   PetscFunctionReturn(0);
442e68589bSMark F. Adams }
452e68589bSMark F. Adams 
462e68589bSMark F. Adams EXTERN_C_BEGIN
472e68589bSMark F. Adams #undef __FUNCT__
482e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG"
492e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n)
502e68589bSMark F. Adams {
512e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
522e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
53c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
542e68589bSMark F. Adams 
552e68589bSMark F. Adams   PetscFunctionBegin;
56c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
57c8b0795cSMark F. Adams   PetscFunctionReturn(0);
58c8b0795cSMark F. Adams }
59c8b0795cSMark F. Adams EXTERN_C_END
60c8b0795cSMark F. Adams 
61c8b0795cSMark F. Adams #undef __FUNCT__
62c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph"
63c8b0795cSMark F. Adams /*@
64c8b0795cSMark F. Adams    PCGAMGSetSymGraph -
65c8b0795cSMark F. Adams 
66c8b0795cSMark F. Adams    Not Collective on PC
67c8b0795cSMark F. Adams 
68c8b0795cSMark F. Adams    Input Parameters:
69c8b0795cSMark F. Adams .  pc - the preconditioner context
70c8b0795cSMark F. Adams 
71c8b0795cSMark F. Adams    Options Database Key:
72c8b0795cSMark F. Adams .  -pc_gamg_sym_graph
73c8b0795cSMark F. Adams 
74c8b0795cSMark F. Adams    Level: intermediate
75c8b0795cSMark F. Adams 
76c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
77c8b0795cSMark F. Adams 
78c8b0795cSMark F. Adams .seealso: ()
79c8b0795cSMark F. Adams @*/
80c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
81c8b0795cSMark F. Adams {
82c8b0795cSMark F. Adams   PetscErrorCode ierr;
83c8b0795cSMark F. Adams 
84c8b0795cSMark F. Adams   PetscFunctionBegin;
85c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
87c8b0795cSMark F. Adams   PetscFunctionReturn(0);
88c8b0795cSMark F. Adams }
89c8b0795cSMark F. Adams 
90c8b0795cSMark F. Adams EXTERN_C_BEGIN
91c8b0795cSMark F. Adams #undef __FUNCT__
92c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph_GAMG"
93c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n)
94c8b0795cSMark F. Adams {
95c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
96c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
97c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
98c8b0795cSMark F. Adams 
99c8b0795cSMark F. Adams   PetscFunctionBegin;
100c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
1012e68589bSMark F. Adams   PetscFunctionReturn(0);
1022e68589bSMark F. Adams }
1032e68589bSMark F. Adams EXTERN_C_END
1042e68589bSMark F. Adams 
105ef4ad70eSMark F. Adams #undef __FUNCT__
106ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph"
107ef4ad70eSMark F. Adams /*@
108ef4ad70eSMark F. Adams    PCGAMGSetSquareGraph -
109ef4ad70eSMark F. Adams 
110ef4ad70eSMark F. Adams    Not Collective on PC
111ef4ad70eSMark F. Adams 
112ef4ad70eSMark F. Adams    Input Parameters:
113ef4ad70eSMark F. Adams .  pc - the preconditioner context
114ef4ad70eSMark F. Adams 
115ef4ad70eSMark F. Adams    Options Database Key:
116ef4ad70eSMark F. Adams .  -pc_gamg_square_graph
117ef4ad70eSMark F. Adams 
118ef4ad70eSMark F. Adams    Level: intermediate
119ef4ad70eSMark F. Adams 
120ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
121ef4ad70eSMark F. Adams 
122ef4ad70eSMark F. Adams .seealso: ()
123ef4ad70eSMark F. Adams @*/
124ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n)
125ef4ad70eSMark F. Adams {
126ef4ad70eSMark F. Adams   PetscErrorCode ierr;
127ef4ad70eSMark F. Adams 
128ef4ad70eSMark F. Adams   PetscFunctionBegin;
129ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
130ef4ad70eSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
131ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
132ef4ad70eSMark F. Adams }
133ef4ad70eSMark F. Adams 
134ef4ad70eSMark F. Adams EXTERN_C_BEGIN
135ef4ad70eSMark F. Adams #undef __FUNCT__
136ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG"
137ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n)
138ef4ad70eSMark F. Adams {
139ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
140ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
141ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
142ef4ad70eSMark F. Adams 
143ef4ad70eSMark F. Adams   PetscFunctionBegin;
144ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
145ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
146ef4ad70eSMark F. Adams }
147ef4ad70eSMark F. Adams EXTERN_C_END
148ef4ad70eSMark F. Adams 
1492e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1502e68589bSMark F. Adams /*
1512e68589bSMark F. Adams    PCSetFromOptions_GAMG_AGG
1522e68589bSMark F. Adams 
1532e68589bSMark F. Adams   Input Parameter:
1542e68589bSMark F. Adams    . pc -
1552e68589bSMark F. Adams */
1562e68589bSMark F. Adams #undef __FUNCT__
1572e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
1582e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc)
1592e68589bSMark F. Adams {
1602e68589bSMark F. Adams   PetscErrorCode ierr;
1612e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1622e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
163c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1642e68589bSMark F. Adams   PetscBool      flag;
1652e68589bSMark F. Adams 
1662e68589bSMark F. Adams   PetscFunctionBegin;
1672e68589bSMark F. Adams   /* call base class */
1682e68589bSMark F. Adams   ierr = PCSetFromOptions_GAMG(pc);CHKERRQ(ierr);
1692e68589bSMark F. Adams 
1702e68589bSMark F. Adams   ierr = PetscOptionsHead("GAMG-AGG options");CHKERRQ(ierr);
1712e68589bSMark F. Adams   {
1722e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
173c8b0795cSMark F. Adams     pc_gamg_agg->nsmooths = 0;
174*2fa5cd67SKarl Rupp 
1752e68589bSMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
1762e68589bSMark F. Adams                            "smoothing steps for smoothed aggregation, usually 1 (0)",
1772e68589bSMark F. Adams                            "PCGAMGSetNSmooths",
178c8b0795cSMark F. Adams                            pc_gamg_agg->nsmooths,
179c8b0795cSMark F. Adams                            &pc_gamg_agg->nsmooths,
180806fa848SBarry Smith                            &flag);CHKERRQ(ierr);
181c8b0795cSMark F. Adams 
182c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
183c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
184*2fa5cd67SKarl Rupp 
185c8b0795cSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_sym_graph",
186581a99e3SJed Brown                             "Set for asymmetric matrices",
187c8b0795cSMark F. Adams                             "PCGAMGSetSymGraph",
188c8b0795cSMark F. Adams                             pc_gamg_agg->sym_graph,
189c8b0795cSMark F. Adams                             &pc_gamg_agg->sym_graph,
190806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
191ef4ad70eSMark F. Adams 
192ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
193ef4ad70eSMark F. Adams     pc_gamg_agg->square_graph = PETSC_TRUE;
194*2fa5cd67SKarl Rupp 
195ef4ad70eSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_square_graph",
1960cbbd2e1SMark F. Adams                             "For faster coarsening and lower coarse grid complexity",
197ef4ad70eSMark F. Adams                             "PCGAMGSetSquareGraph",
198ef4ad70eSMark F. Adams                             pc_gamg_agg->square_graph,
199ef4ad70eSMark F. Adams                             &pc_gamg_agg->square_graph,
200806fa848SBarry Smith                             &flag);CHKERRQ(ierr);
2012e68589bSMark F. Adams   }
2022e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
2032e68589bSMark F. Adams   PetscFunctionReturn(0);
2042e68589bSMark F. Adams }
2052e68589bSMark F. Adams 
2062e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2072e68589bSMark F. Adams /*
2082e68589bSMark F. Adams    PCDestroy_AGG
2092e68589bSMark F. Adams 
2102e68589bSMark F. Adams   Input Parameter:
2112e68589bSMark F. Adams    . pc -
2122e68589bSMark F. Adams */
2132e68589bSMark F. Adams #undef __FUNCT__
2142e68589bSMark F. Adams #define __FUNCT__ "PCDestroy_AGG"
2152e68589bSMark F. Adams PetscErrorCode PCDestroy_AGG(PC pc)
2162e68589bSMark F. Adams {
2172e68589bSMark F. Adams   PetscErrorCode ierr;
2182e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
2192e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
220c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
2212e68589bSMark F. Adams 
2222e68589bSMark F. Adams   PetscFunctionBegin;
223c8b0795cSMark F. Adams   if (pc_gamg_agg) {
224c8b0795cSMark F. Adams     ierr        = PetscFree(pc_gamg_agg);CHKERRQ(ierr);
225c8b0795cSMark F. Adams     pc_gamg_agg = 0;
2262e68589bSMark F. Adams   }
2272e68589bSMark F. Adams 
2282e68589bSMark F. Adams   /* call base class */
2292e68589bSMark F. Adams   ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
2302e68589bSMark F. Adams   PetscFunctionReturn(0);
2312e68589bSMark F. Adams }
2322e68589bSMark F. Adams 
2332e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2342e68589bSMark F. Adams /*
2352e68589bSMark F. Adams    PCSetCoordinates_AGG
236302f38e8SMark F. Adams      - collective
2372e68589bSMark F. Adams 
2382e68589bSMark F. Adams    Input Parameter:
2392e68589bSMark F. Adams    . pc - the preconditioner context
240a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
241302f38e8SMark F. Adams    . a_nloc - number of vertices local
242302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2432e68589bSMark F. Adams */
2442e68589bSMark F. Adams EXTERN_C_BEGIN
2452e68589bSMark F. Adams #undef __FUNCT__
2462e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
247302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2482e68589bSMark F. Adams {
2492e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2502e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2512e68589bSMark F. Adams   PetscErrorCode ierr;
25269344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
253a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
254a2f3521dSMark F. Adams   /* MPI_Comm     wcomm = ((PetscObject)pc)->comm; */
2552e68589bSMark F. Adams 
2562e68589bSMark F. Adams   PetscFunctionBegin;
257a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
258a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
259302f38e8SMark F. Adams   nloc = a_nloc;
2602e68589bSMark F. Adams 
2612e68589bSMark F. Adams   /* SA: null space vectors */
26269344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
26369344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
264a2f3521dSMark F. Adams   else if (coords) {
265c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
26669344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
26769344418SMark F. Adams     if (ndm != ndf) {
268c666adbfSMark 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);
269a2f3521dSMark F. Adams     }
270806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
27169344418SMark F. Adams   pc_gamg->data_cell_rows = ndatarows = ndf;     assert(pc_gamg->data_cell_cols>0);
272*2fa5cd67SKarl Rupp 
273c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2742e68589bSMark F. Adams 
2752e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2762e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2772e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
278302f38e8SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
2792e68589bSMark F. Adams   }
2802e68589bSMark F. Adams   /* copy data in - column oriented */
2812e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
28269344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
28369344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
284c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2852e68589bSMark F. Adams     else {
28669344418SMark F. Adams       /* translational modes */
287*2fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
288*2fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
28969344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2902e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
291*2fa5cd67SKarl Rupp         }
292*2fa5cd67SKarl Rupp       }
29369344418SMark F. Adams 
29469344418SMark F. Adams       /* rotational modes */
2952e68589bSMark F. Adams       if (coords) {
29669344418SMark F. Adams         if (ndm == 2) {
2972e68589bSMark F. Adams           data   += 2*M;
2982e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2992e68589bSMark F. Adams           data[1] =  coords[2*kk];
300806fa848SBarry Smith         } else {
3012e68589bSMark F. Adams           data   += 3*M;
3022e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
3032e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
3042e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
3052e68589bSMark F. Adams         }
3062e68589bSMark F. Adams       }
3072e68589bSMark F. Adams     }
3082e68589bSMark F. Adams   }
3092e68589bSMark F. Adams 
3102e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
3112e68589bSMark F. Adams   PetscFunctionReturn(0);
3122e68589bSMark F. Adams }
3132e68589bSMark F. Adams EXTERN_C_END
3142e68589bSMark F. Adams 
315b43b03e9SMark F. Adams typedef PetscInt NState;
316b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
317b43b03e9SMark F. Adams static const NState DELETED =-1;
318b43b03e9SMark F. Adams static const NState REMOVED =-3;
319b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
320b43b03e9SMark F. Adams 
321c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
322c8b0795cSMark F. Adams /*
323b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
324b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
325c8b0795cSMark F. Adams 
326c8b0795cSMark F. Adams    Input Parameter:
327c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
328c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
329c8b0795cSMark F. Adams    Input/Output Parameter:
3300cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
331c8b0795cSMark F. Adams */
332c8b0795cSMark F. Adams #undef __FUNCT__
333c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3340cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
3350cbbd2e1SMark F. Adams                                  const Mat Gmat_1,  /* base graph */
3360cbbd2e1SMark F. Adams                                  /* const IS selected_2, [nselected local] selected vertices */
3371147fc2aSKarl Rupp                                  PetscCoarsenData *aggs_2)  /* [nselected local] global ID of aggregate */
338c8b0795cSMark F. Adams {
339c8b0795cSMark F. Adams   PetscErrorCode ierr;
340c8b0795cSMark F. Adams   PetscBool      isMPI;
341c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
342c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat_2)->comm;
343c5df96a5SBarry Smith   PetscMPIInt    rank,size;
3440cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
345c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
346c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3470cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3480cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
349c8b0795cSMark F. Adams   NState         *lid_state;
3500cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
351c8b0795cSMark F. Adams 
352c8b0795cSMark F. Adams   PetscFunctionBegin;
353c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
354c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
355c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
356c8b0795cSMark F. Adams 
3570cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
358c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
359c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
360c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
361c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
362c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
363c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
364c8b0795cSMark F. Adams   }
365c8b0795cSMark F. Adams 
366c8b0795cSMark F. Adams   /* get submatrices */
367251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
368c8b0795cSMark F. Adams   if (isMPI) {
369c8b0795cSMark F. Adams     /* grab matrix objects */
370c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
371c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
372c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
373c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
374c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
375c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
376c8b0795cSMark F. Adams 
377c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
378c8b0795cSMark F. Adams     matB_1->compressedrow.check = PETSC_TRUE;
379*2fa5cd67SKarl Rupp 
380806fa848SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
381c8b0795cSMark F. Adams 
382c8b0795cSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &lid_cprowID_1);CHKERRQ(ierr);
3830cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
384c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
385c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
386c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
387c8b0795cSMark F. Adams     }
388806fa848SBarry Smith   } else {
389c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
390c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3910cbbd2e1SMark F. Adams     lid_cprowID_1 = PETSC_NULL;
392c8b0795cSMark F. Adams   }
393c8b0795cSMark F. Adams   assert(matA_1 && !matA_1->compressedrow.use);
394c8b0795cSMark F. Adams   assert(matB_1==0 || matB_1->compressedrow.use);
395c8b0795cSMark F. Adams   assert(matA_2 && !matA_2->compressedrow.use);
396c8b0795cSMark F. Adams   assert(matB_2==0 || matB_2->compressedrow.use);
397c8b0795cSMark F. Adams 
398c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
399c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(NState), &lid_state);CHKERRQ(ierr);
4000cbbd2e1SMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscScalar), &lid_parent_gid);CHKERRQ(ierr);
401c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
4020cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
403c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
404c8b0795cSMark F. Adams   }
4050cbbd2e1SMark F. Adams 
4060cbbd2e1SMark F. Adams   /* set lid_state */
4070cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
40841b27cdeSMark F. Adams     PetscCDPos pos;
409e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
410e78576d6SMark F. Adams     if (pos) {
411e78576d6SMark F. Adams       PetscInt gid1;
412*2fa5cd67SKarl Rupp 
413ffc955d6SMark F. Adams       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); assert(gid1==lid+my0);
414*2fa5cd67SKarl Rupp 
4150cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
416b43b03e9SMark F. Adams     }
417b43b03e9SMark F. Adams   }
4180cbbd2e1SMark F. Adams 
4190cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
420c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
421c8b0795cSMark F. Adams     NState state = lid_state[lid];
422c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
42341b27cdeSMark F. Adams       PetscCDPos pos;
424e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
425e78576d6SMark F. Adams       while (pos) {
426e78576d6SMark F. Adams         PetscInt gid1;
427ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
428e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
429e78576d6SMark F. Adams 
430*2fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
431c8b0795cSMark F. Adams       }
4320cbbd2e1SMark F. Adams     }
4330cbbd2e1SMark F. Adams   }
4340cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
435c8b0795cSMark F. Adams   if (isMPI) {
436c8b0795cSMark F. Adams     Vec tempVec;
437c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
438c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
439c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
440c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
441c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
442c8b0795cSMark F. Adams     }
443c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
444c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
445806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
446806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
447c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
448c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
449806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
450806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
451c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4520cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4530cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4540cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4550cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4560cbbd2e1SMark F. Adams     }
4570cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4580cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4590cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
460806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
461806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4620cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4630cbbd2e1SMark F. Adams 
464c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
465c8b0795cSMark F. Adams   } /* ismpi */
466c8b0795cSMark F. Adams 
467c8b0795cSMark F. Adams   /* doit */
468c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
469c8b0795cSMark F. Adams     NState state = lid_state[lid];
4700cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4710cbbd2e1SMark F. Adams       /* steal locals */
472c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
473c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
474c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4750cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
476c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4770cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4780cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4790cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4800cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
48141b27cdeSMark F. Adams             PetscCDPos pos,last=PETSC_NULL;
482c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
483e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
484e78576d6SMark F. Adams             while (pos) {
485e78576d6SMark F. Adams               PetscInt gid;
486ffc955d6SMark F. Adams               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
4870cbbd2e1SMark F. Adams               if (gid == gidj) {
4880cbbd2e1SMark F. Adams                 assert(last);
48941b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
49041b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4910cbbd2e1SMark F. Adams                 hav  = 1;
492c8b0795cSMark F. Adams                 break;
493806fa848SBarry Smith               } else last = pos;
494e78576d6SMark F. Adams 
495e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
496c8b0795cSMark F. Adams             }
497c8b0795cSMark F. Adams             if (hav!=1) {
498c666adbfSMark F. Adams               if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
499c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
500c8b0795cSMark F. Adams             }
501806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
502*2fa5cd67SKarl 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.");
50341b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
504c8b0795cSMark F. Adams           }
505c8b0795cSMark F. Adams         }
5060cbbd2e1SMark F. Adams       } /* local neighbors */
507806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
5080cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
509c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
510c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
511c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
512c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
513c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
514c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
515c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
516c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
5170cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
5180cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
5190cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
52041b27cdeSMark F. Adams               PetscCDPos pos,last=PETSC_NULL;
5210cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
522e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
523e78576d6SMark F. Adams               while (pos) {
524e78576d6SMark F. Adams                 PetscInt gid;
525ffc955d6SMark F. Adams                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
5260cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
5270cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
5280cbbd2e1SMark F. Adams                   assert(last);
52941b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
5300cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
5310cbbd2e1SMark F. Adams                   hav = 1;
532c8b0795cSMark F. Adams                   break;
533806fa848SBarry Smith                 } else last = pos;
534e78576d6SMark F. Adams 
535e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
536c8b0795cSMark F. Adams               }
537c8b0795cSMark F. Adams               if (hav!=1) {
538c666adbfSMark F. Adams                 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
539c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
540c8b0795cSMark F. Adams               }
541806fa848SBarry Smith             } else {
5420cbbd2e1SMark F. Adams               /* ghosts remove this later */
5430cbbd2e1SMark F. Adams             }
544c8b0795cSMark F. Adams           }
545c8b0795cSMark F. Adams         }
546c8b0795cSMark F. Adams       }
547c8b0795cSMark F. Adams     } /* selected/deleted */
548c8b0795cSMark F. Adams   } /* node loop */
549c8b0795cSMark F. Adams 
550c8b0795cSMark F. Adams   if (isMPI) {
5510cbbd2e1SMark F. Adams     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
5520cbbd2e1SMark F. Adams     Vec           tempVec,ghostgids2,ghostparents2;
5530cbbd2e1SMark F. Adams     PetscInt      cpid,nghost_2;
5540cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
555c8b0795cSMark F. Adams 
5560cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
557c8b0795cSMark F. Adams     ierr = MatGetVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5580cbbd2e1SMark F. Adams 
5590cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
560c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5610cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
562c8b0795cSMark F. Adams     }
563c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
564c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5650cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
566806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
567806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5680cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5690cbbd2e1SMark F. Adams 
5700cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5710cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5720cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5730cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5740cbbd2e1SMark F. Adams     }
5750cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5760cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5770cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
578806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
579806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5800cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5810cbbd2e1SMark F. Adams 
582c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
583c8b0795cSMark F. Adams 
5840cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
585f10ff945SMark F. Adams     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5860cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5870cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5880cbbd2e1SMark F. Adams       if (state==DELETED) {
5890cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5900cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5910cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5920cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5930cbbd2e1SMark F. Adams           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5940cbbd2e1SMark F. Adams         }
5950cbbd2e1SMark F. Adams       }
5960cbbd2e1SMark F. Adams     }
597c8b0795cSMark F. Adams 
5980cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
599c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
600c8b0795cSMark F. Adams       NState state = lid_state[lid];
601c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
60241b27cdeSMark F. Adams         PetscCDPos pos,last=PETSC_NULL;
603c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
604e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
605e78576d6SMark F. Adams         while (pos) {
606e78576d6SMark F. Adams           PetscInt gid;
607ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
608e78576d6SMark F. Adams 
6090cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
6100cbbd2e1SMark F. Adams             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
6110cbbd2e1SMark F. Adams             if (cpid != -1) {
6120cbbd2e1SMark F. Adams               /* a moved ghost - */
6130cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
61441b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
615806fa848SBarry Smith             } else last = pos;
616806fa848SBarry Smith           } else last = pos;
617e78576d6SMark F. Adams 
618e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
619c8b0795cSMark F. Adams         } /* loop over list of deleted */
620c8b0795cSMark F. Adams       } /* selected */
621c8b0795cSMark F. Adams     }
6220cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
623c8b0795cSMark F. Adams 
6240cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
6250cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
6260cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6270cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
6280cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6290cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
63041b27cdeSMark F. Adams         PetscCDPos pos;
6310cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
632e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
633e78576d6SMark F. Adams         while (pos) {
634e78576d6SMark F. Adams           PetscInt gidj;
635ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
636e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
637e78576d6SMark F. Adams 
6380cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
639c8b0795cSMark F. Adams         }
640c8b0795cSMark F. Adams         if (hav != 1) {
641ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
64241b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
643c8b0795cSMark F. Adams         }
644c8b0795cSMark F. Adams       }
645c8b0795cSMark F. Adams     }
646c8b0795cSMark F. Adams 
6470cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6480cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6490cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6500cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
651c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6520cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6530cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6540cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
655c8b0795cSMark F. Adams   }
656c8b0795cSMark F. Adams 
6570cbbd2e1SMark F. Adams   ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr);
658c8b0795cSMark F. Adams   ierr = PetscFree(lid_state);CHKERRQ(ierr);
659c8b0795cSMark F. Adams   PetscFunctionReturn(0);
660c8b0795cSMark F. Adams }
6612e68589bSMark F. Adams 
6622e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6632e68589bSMark F. Adams /*
664a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
665a2f3521dSMark F. Adams       Looks in Mat for near null space.
666a2f3521dSMark F. Adams       Does not work for Stokes
6672e68589bSMark F. Adams 
6682e68589bSMark F. Adams   Input Parameter:
6692e68589bSMark F. Adams    . pc -
670a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6712e68589bSMark F. Adams */
6722e68589bSMark F. Adams #undef __FUNCT__
6732e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
674b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6752e68589bSMark F. Adams {
6762e68589bSMark F. Adams   PetscErrorCode ierr;
677b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
678b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
679b8cd405aSMark F. Adams   MatNullSpace   mnull;
680b8cd405aSMark F. Adams 
6812e68589bSMark F. Adams   PetscFunctionBegin;
682b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
683b8cd405aSMark F. Adams   if (!mnull) {
684a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6859d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
6869d1b15c3SMark F. Adams     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);  assert(MM%bs==0);
687a2f3521dSMark F. Adams     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, PETSC_NULL);CHKERRQ(ierr);
688806fa848SBarry Smith   } else {
689b8cd405aSMark F. Adams     PetscReal *nullvec;
690b8cd405aSMark F. Adams     PetscBool has_const;
691b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
692b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
693b8cd405aSMark F. Adams     ierr = MatGetLocalSize(a_A,&mlocal,PETSC_NULL);CHKERRQ(ierr);
694b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
6958caf3d72SBarry Smith     ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
696b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
697b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
698b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
699b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
700b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
701b8cd405aSMark F. Adams     }
702b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
703b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
704*2fa5cd67SKarl Rupp 
705a2f3521dSMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); /* this does not work for Stokes */
706*2fa5cd67SKarl Rupp 
707b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
708b8cd405aSMark F. Adams   }
7092e68589bSMark F. Adams   PetscFunctionReturn(0);
7102e68589bSMark F. Adams }
7112e68589bSMark F. Adams 
7122e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
7132e68589bSMark F. Adams /*
7142e68589bSMark F. Adams  formProl0
7152e68589bSMark F. Adams 
7162e68589bSMark F. Adams    Input Parameter:
7170cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
7182e68589bSMark F. Adams    . bs - block size
7190cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7200cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7212e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
7222e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7232e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7242e68589bSMark F. Adams   Output Parameter:
7252e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7262e68589bSMark F. Adams    . a_Prol - prolongation operator
7272e68589bSMark F. Adams */
7282e68589bSMark F. Adams #undef __FUNCT__
7292e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7300cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */
7310cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
7320cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
7330cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
7340cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7350cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7360cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7372e68589bSMark F. Adams                                 PetscReal **a_data_out,
7381147fc2aSKarl Rupp                                 Mat a_Prol) /* prolongation operator (output)*/
7392e68589bSMark F. Adams {
7402e68589bSMark F. Adams   PetscErrorCode ierr;
7410cbbd2e1SMark F. Adams   PetscInt       Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7422e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
743c5df96a5SBarry Smith   PetscMPIInt    rank, size;
7442e68589bSMark F. Adams   PetscReal      *out_data;
74541b27cdeSMark F. Adams   PetscCDPos     pos;
7460cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7470cbbd2e1SMark F. Adams 
748797e13b7SMark F. Adams /* #define OUT_AGGS */
749519f805aSKarl Rupp #if defined(OUT_AGGS)
750f7620de1SMatthew G Knepley   static PetscInt llev = 0; char fname[32]; FILE *file = PETSC_NULL; PetscInt pM;
7519057884aSMark F. Adams #endif
7522e68589bSMark F. Adams 
7532e68589bSMark F. Adams   PetscFunctionBegin;
754c5df96a5SBarry Smith   ierr    = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
755c5df96a5SBarry Smith   ierr    = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
7562e68589bSMark F. Adams   ierr    = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
7570cbbd2e1SMark F. Adams   nloc    = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
7580cbbd2e1SMark F. Adams   Iend   /= bs;
7590cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7602e68589bSMark F. Adams 
761f10ff945SMark F. Adams   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7620cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7630cbbd2e1SMark F. Adams     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7642e68589bSMark F. Adams   }
7652e68589bSMark F. Adams 
766519f805aSKarl Rupp #if defined(OUT_AGGS)
767c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
768*2fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7690cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7700cbbd2e1SMark F. Adams #endif
7710cbbd2e1SMark F. Adams 
7720cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7730cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
774e78576d6SMark F. Adams     PetscBool ise;
775e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
776e78576d6SMark F. Adams     if (!ise) nSelected++;
7770cbbd2e1SMark F. Adams   }
7780cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
7790cbbd2e1SMark F. Adams   assert((ii/nSAvec)==my0crs); assert(nSelected==(jj-ii)/nSAvec);
7800cbbd2e1SMark F. Adams 
7812e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7820cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
783*2fa5cd67SKarl Rupp 
7840cbbd2e1SMark F. Adams   ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr);
785*2fa5cd67SKarl Rupp   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=1.e300;
7860cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7872e68589bSMark F. Adams 
7882e68589bSMark F. Adams   /* find points and set prolongation */
789c8b0795cSMark F. Adams   minsz = 100;
7902e68589bSMark F. Adams   ndone = 0;
7910cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
792e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
793e78576d6SMark F. Adams     if (jj > 0) {
7940cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7950cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7960cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7972e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7982e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7992e68589bSMark F. Adams       PetscInt       *fids;
80065d7b583SSatish Balay       PetscReal      *data;
8010cbbd2e1SMark F. Adams       /* count agg */
8020cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
8030cbbd2e1SMark F. Adams 
8040cbbd2e1SMark F. Adams       /* get block */
8052e68589bSMark F. Adams       ierr = PetscMalloc((Mdata*N)*sizeof(PetscScalar), &qqc);CHKERRQ(ierr);
8062e68589bSMark F. Adams       ierr = PetscMalloc((M*N)*sizeof(PetscScalar), &qqr);CHKERRQ(ierr);
8072e68589bSMark F. Adams       ierr = PetscMalloc(N*sizeof(PetscScalar), &TAU);CHKERRQ(ierr);
8082e68589bSMark F. Adams       ierr = PetscMalloc(LWORK*sizeof(PetscScalar), &WORK);CHKERRQ(ierr);
8092e68589bSMark F. Adams       ierr = PetscMalloc(M*sizeof(PetscInt), &fids);CHKERRQ(ierr);
8102e68589bSMark F. Adams 
8112e68589bSMark F. Adams       aggID = 0;
812e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
813e78576d6SMark F. Adams       while (pos) {
814e78576d6SMark F. Adams         PetscInt gid1;
815ffc955d6SMark F. Adams         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
816e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
817e78576d6SMark F. Adams 
8180cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
8190cbbd2e1SMark F. Adams         else {
8200cbbd2e1SMark F. Adams           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
8210cbbd2e1SMark F. Adams           assert(flid>=0);
8220cbbd2e1SMark F. Adams         }
8232e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
82465d7b583SSatish Balay         data = &data_in[flid*bs];
8252e68589bSMark F. Adams         for (kk = ii = 0; ii < bs; ii++) {
8262e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
8270cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8280cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8292e68589bSMark F. Adams           }
8302e68589bSMark F. Adams         }
831519f805aSKarl Rupp #if defined(OUT_AGGS)
832b2a4f308SMark F. Adams         if (llev==1) {
8339057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8340cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
835c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
836f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8370cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
838b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
839b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8409057884aSMark F. Adams         }
8419057884aSMark F. Adams #endif
8422e68589bSMark F. Adams         /* set fine IDs */
8432e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8442e68589bSMark F. Adams 
8452e68589bSMark F. Adams         aggID++;
8460cbbd2e1SMark F. Adams       }
8472e68589bSMark F. Adams 
8482e68589bSMark F. Adams       /* pad with zeros */
8492e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8502e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8512e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8522e68589bSMark F. Adams         }
8532e68589bSMark F. Adams       }
8542e68589bSMark F. Adams 
8552e68589bSMark F. Adams       ndone += aggID;
8562e68589bSMark F. Adams       /* QR */
85784df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8582e68589bSMark F. Adams       LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO);
85984df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
860d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8612e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8622e68589bSMark F. Adams       {
8632e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8642e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8652e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
8660cbbd2e1SMark F. Adams             assert(data[jj*out_data_stride + ii] == 1.e300);
8670cbbd2e1SMark F. Adams             if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8680cbbd2e1SMark F. Adams             else data[jj*out_data_stride + ii] = 0.;
8692e68589bSMark F. Adams           }
8702e68589bSMark F. Adams         }
8712e68589bSMark F. Adams       }
8722e68589bSMark F. Adams 
8732e68589bSMark F. Adams       /* get Q - row oriented */
8742e68589bSMark F. Adams       LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO);
875c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8762e68589bSMark F. Adams 
8772e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8782e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8792e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8802e68589bSMark F. Adams         }
8812e68589bSMark F. Adams       }
8822e68589bSMark F. Adams 
8832e68589bSMark F. Adams       /* add diagonal block of P0 */
884c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
885c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
886c8b0795cSMark F. Adams       }
8872e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8882e68589bSMark F. Adams 
8892e68589bSMark F. Adams       ierr = PetscFree(qqc);CHKERRQ(ierr);
8902e68589bSMark F. Adams       ierr = PetscFree(qqr);CHKERRQ(ierr);
8912e68589bSMark F. Adams       ierr = PetscFree(TAU);CHKERRQ(ierr);
8922e68589bSMark F. Adams       ierr = PetscFree(WORK);CHKERRQ(ierr);
8932e68589bSMark F. Adams       ierr = PetscFree(fids);CHKERRQ(ierr);
894b43b03e9SMark F. Adams       clid++;
8950cbbd2e1SMark F. Adams     } /* coarse agg */
8960cbbd2e1SMark F. Adams   } /* for all fine nodes */
8970cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8980cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8992e68589bSMark F. Adams 
900c8b0795cSMark F. Adams /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm); */
9012e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
902c8b0795cSMark F. Adams /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm); */
903c5df96a5SBarry Smith /* PetscPrintf(wcomm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */
9042e68589bSMark F. Adams 
905519f805aSKarl Rupp #if defined(OUT_AGGS)
906b2a4f308SMark F. Adams   if (llev==1) fclose(file);
9079057884aSMark F. Adams #endif
9080cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
9092e68589bSMark F. Adams   PetscFunctionReturn(0);
9102e68589bSMark F. Adams }
9112e68589bSMark F. Adams 
9122e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
9132e68589bSMark F. Adams /*
914c8b0795cSMark F. Adams    PCGAMGgraph_AGG
9152e68589bSMark F. Adams 
9162e68589bSMark F. Adams   Input Parameter:
9172e68589bSMark F. Adams    . pc - this
9182e68589bSMark F. Adams    . Amat - matrix on this fine level
9192e68589bSMark F. Adams   Output Parameter:
920c8b0795cSMark F. Adams    . a_Gmat -
9212e68589bSMark F. Adams */
9222e68589bSMark F. Adams #undef __FUNCT__
923c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG"
924*2fa5cd67SKarl Rupp PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat)
925c8b0795cSMark F. Adams {
926c8b0795cSMark F. Adams   PetscErrorCode            ierr;
927c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
928c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
929c8b0795cSMark F. Adams   const PetscInt            verbose      = pc_gamg->verbose;
930c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
931c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
932c5df96a5SBarry Smith   PetscMPIInt               rank,size;
933e0940f08SMark F. Adams   Mat                       Gmat;
934c8b0795cSMark F. Adams   MPI_Comm                  wcomm = ((PetscObject)Amat)->comm;
93567744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
936c8b0795cSMark F. Adams 
937c8b0795cSMark F. Adams   PetscFunctionBegin;
9380cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9390cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9400cbbd2e1SMark F. Adams #endif
941c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
942c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
943c8b0795cSMark F. Adams 
94467744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
945c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9460cbbd2e1SMark F. Adams 
9472d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
9482d7fac45SMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
949c8b0795cSMark F. Adams 
950e0940f08SMark F. Adams   *a_Gmat = Gmat;
951c8b0795cSMark F. Adams 
9520cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9530cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9540cbbd2e1SMark F. Adams #endif
955c8b0795cSMark F. Adams   PetscFunctionReturn(0);
956c8b0795cSMark F. Adams }
957c8b0795cSMark F. Adams 
958c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
959c8b0795cSMark F. Adams /*
960b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
961c8b0795cSMark F. Adams 
962c8b0795cSMark F. Adams   Input Parameter:
963e0940f08SMark F. Adams    . a_pc - this
964e0940f08SMark F. Adams   Input/Output Parameter:
9650cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
966c8b0795cSMark F. Adams   Output Parameter:
9670cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
968c8b0795cSMark F. Adams */
969c8b0795cSMark F. Adams #undef __FUNCT__
970b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
971*2fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
972c8b0795cSMark F. Adams {
973c8b0795cSMark F. Adams   PetscErrorCode ierr;
974e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
975c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
976c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9770cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9789f3f12c8SMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
9790cbbd2e1SMark F. Adams   IS             perm;
980c8b0795cSMark F. Adams   PetscInt       Ii,nloc,bs,n,m;
981c8b0795cSMark F. Adams   PetscInt       *permute;
982c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
983b43b03e9SMark F. Adams   MatCoarsen     crs;
984e0940f08SMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
985c5df96a5SBarry Smith   PetscMPIInt    rank,size;
986c8b0795cSMark F. Adams 
987c8b0795cSMark F. Adams   PetscFunctionBegin;
9880cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9890cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9900cbbd2e1SMark F. Adams #endif
991c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
992c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
993e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
994e0940f08SMark F. Adams   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); assert(bs==1);
995c8b0795cSMark F. Adams   nloc = n/bs;
996c8b0795cSMark F. Adams 
997e0940f08SMark F. Adams   if (pc_gamg_agg->square_graph) {
998*2fa5cd67SKarl Rupp     if (verbose > 1) PetscPrintf(wcomm,"[%d]%s square graph\n",rank,__FUNCT__);
999878e152fSMark F. Adams     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
1000806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
10019f3f12c8SMark F. Adams     if (verbose > 2) {
1002c5df96a5SBarry Smith       ierr = PetscPrintf(wcomm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr);
1003f10ff945SMark F. Adams       /* check for symetry */
1004f10ff945SMark F. Adams       if (verbose > 4) {
1005f10ff945SMark F. Adams 
1006f10ff945SMark F. Adams       }
10079f3f12c8SMark F. Adams     }
1008806fa848SBarry Smith   } else Gmat2 = Gmat1;
1009c8b0795cSMark F. Adams 
1010c8b0795cSMark F. Adams   /* get MIS aggs */
1011c8b0795cSMark F. Adams   /* randomize */
1012c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
1013c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
1014c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1015c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
1016c8b0795cSMark F. Adams     permute[Ii]   = Ii;
1017c8b0795cSMark F. Adams   }
1018c8b0795cSMark F. Adams   srand(1); /* make deterministic */
1019c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1020c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
1021c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1022c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1023c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1024c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1025c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1026c8b0795cSMark F. Adams     }
1027c8b0795cSMark F. Adams   }
1028c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
1029c8b0795cSMark F. Adams 
1030*2fa5cd67SKarl Rupp   if (verbose > 1) PetscPrintf(wcomm,"[%d]%s coarsen graph\n",rank,__FUNCT__);
10319f3f12c8SMark F. Adams 
1032806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
10330cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10340cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1035b43b03e9SMark F. Adams #endif
1036b43b03e9SMark F. Adams   ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
10379057884aSMark F. Adams   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
10389057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1039b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1040b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1041b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
1042b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1043b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
10440cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1045b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1046b43b03e9SMark F. Adams 
1047c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1048c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10500cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1051b43b03e9SMark F. Adams #endif
1052*2fa5cd67SKarl Rupp   if (verbose > 2) PetscPrintf(wcomm,"[%d]%s coarsen graph done\n",rank,__FUNCT__);
10539f3f12c8SMark F. Adams 
1054c8b0795cSMark F. Adams   /* smooth aggs */
1055e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10560cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10570cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1058c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1059e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
106041b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10610cbbd2e1SMark F. Adams     if (mat) SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1062806fa848SBarry Smith   } else {
10630cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1064f10ff945SMark F. Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
106541b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10660cbbd2e1SMark F. Adams     if (mat) {
10670cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10680cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10690cbbd2e1SMark F. Adams     }
10700cbbd2e1SMark F. Adams   }
10710cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
10720cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
10730cbbd2e1SMark F. Adams #endif
1074*2fa5cd67SKarl Rupp   if (verbose > 2) PetscPrintf(wcomm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__);
1075c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1076c8b0795cSMark F. Adams }
1077c8b0795cSMark F. Adams 
1078c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1079c8b0795cSMark F. Adams /*
10800cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1081c8b0795cSMark F. Adams 
1082c8b0795cSMark F. Adams  Input Parameter:
1083c8b0795cSMark F. Adams  . pc - this
1084c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1085c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10860cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1087c8b0795cSMark F. Adams  Output Parameter:
1088c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1089c8b0795cSMark F. Adams  */
1090c8b0795cSMark F. Adams #undef __FUNCT__
10910cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
1092*2fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10932e68589bSMark F. Adams {
10942e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10952e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
10962e68589bSMark F. Adams   const PetscInt verbose   = pc_gamg->verbose;
1097c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10982e68589bSMark F. Adams   PetscErrorCode ierr;
1099c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1100c8b0795cSMark F. Adams   Mat            Prol;
1101c5df96a5SBarry Smith   PetscMPIInt    rank, size;
11022e68589bSMark F. Adams   MPI_Comm       wcomm  = ((PetscObject)Amat)->comm;
11030cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1104c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1105c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
11062e68589bSMark F. Adams 
11072e68589bSMark F. Adams   PetscFunctionBegin;
11080cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11090cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11100cbbd2e1SMark F. Adams #endif
1111c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
1112c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
11132e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1114c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
1115c8b0795cSMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
11162e68589bSMark F. Adams 
11172e68589bSMark F. Adams   /* get 'nLocalSelected' */
11180cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1119e78576d6SMark F. Adams     PetscBool ise;
11200cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1121e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1122e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
11232e68589bSMark F. Adams   }
11242e68589bSMark F. Adams 
11252e68589bSMark F. Adams   /* create prolongator, create P matrix */
1126a2f3521dSMark F. Adams   ierr = MatCreate(wcomm, &Prol);CHKERRQ(ierr);
1127806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1128a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1129a2f3521dSMark F. Adams   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
1130a2f3521dSMark F. Adams   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, PETSC_NULL);CHKERRQ(ierr);
1131a2f3521dSMark F. Adams   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, PETSC_NULL,data_cols, PETSC_NULL);CHKERRQ(ierr);
1132a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1133a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
1134a2f3521dSMark F. Adams   /* data_cols, PETSC_NULL, data_cols, PETSC_NULL, */
1135a2f3521dSMark F. Adams   /* &Prol); */
11362e68589bSMark F. Adams 
11372e68589bSMark F. Adams   /* can get all points "removed" */
1138c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1139c8b0795cSMark F. Adams   if (ii==0) {
1140*2fa5cd67SKarl Rupp     if (verbose > 0) PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__);
11412e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
11422e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
11432e68589bSMark F. Adams     PetscFunctionReturn(0);
11442e68589bSMark F. Adams   }
1145*2fa5cd67SKarl Rupp   if (verbose > 0) PetscPrintf(wcomm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs);
1146c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
11470cbbd2e1SMark F. Adams 
11480cbbd2e1SMark F. Adams   assert((kk-myCrs0)%col_bs==0);
1149c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
11500cbbd2e1SMark F. Adams   assert((kk/col_bs-myCrs0)==nLocalSelected);
11512e68589bSMark F. Adams 
11522e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11530cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11540cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11552e68589bSMark F. Adams #endif
1156c5df96a5SBarry Smith   if (size > 1) { /*  */
11572e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
11582e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr);
11592e68589bSMark F. Adams     for (jj = 0; jj < data_cols; jj++) {
1160c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1161a2f3521dSMark F. Adams         PetscInt        ii,stride;
1162c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1163*2fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
1164*2fa5cd67SKarl Rupp 
1165806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1166a2f3521dSMark F. Adams 
11672e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1168a2f3521dSMark F. Adams           ierr    = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr);
1169a2f3521dSMark F. Adams           nbnodes = bs*stride;
11702e68589bSMark F. Adams         }
1171a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
1172*2fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11732e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11742e68589bSMark F. Adams       }
11752e68589bSMark F. Adams     }
11762e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1177806fa848SBarry Smith   } else {
1178c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1179c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11802e68589bSMark F. Adams   }
11812e68589bSMark F. Adams 
11822e68589bSMark F. Adams   /* get P0 */
1183c5df96a5SBarry Smith   if (size > 1) {
11842e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1185a2f3521dSMark F. Adams     PetscInt  stride;
11862e68589bSMark F. Adams 
11872e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr);
11882e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1189806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1190a2f3521dSMark F. Adams     ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1191a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11922e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1193a2f3521dSMark F. Adams 
1194a2f3521dSMark F. Adams     assert(stride==nbnodes/bs);
11952e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1196806fa848SBarry Smith   } else {
11972e68589bSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
11982e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11992e68589bSMark F. Adams   }
12000cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12010cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
12020cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
12032e68589bSMark F. Adams #endif
1204c8b0795cSMark F. Adams   {
1205ffc955d6SMark F. Adams     PetscReal *data_out = PETSC_NULL;
1206806fa848SBarry Smith     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1207c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
1208*2fa5cd67SKarl Rupp 
1209c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
1210c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1211c8b0795cSMark F. Adams     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1212c8b0795cSMark F. Adams   }
12130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12140cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1215c8b0795cSMark F. Adams #endif
1216c5df96a5SBarry Smith   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
12172e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
12182e68589bSMark F. Adams 
1219c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
12200cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12210cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
12220cbbd2e1SMark F. Adams #endif
1223c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1224c8b0795cSMark F. Adams }
1225c8b0795cSMark F. Adams 
1226c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1227c8b0795cSMark F. Adams /*
12280cbbd2e1SMark F. Adams    PCGAMGOptprol_AGG
1229c8b0795cSMark F. Adams 
1230c8b0795cSMark F. Adams   Input Parameter:
1231c8b0795cSMark F. Adams    . pc - this
1232c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1233c8b0795cSMark F. Adams  In/Output Parameter:
1234c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1235c8b0795cSMark F. Adams */
1236c8b0795cSMark F. Adams #undef __FUNCT__
12370cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG"
1238*2fa5cd67SKarl Rupp PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P)
1239c8b0795cSMark F. Adams {
1240c8b0795cSMark F. Adams   PetscErrorCode ierr;
1241c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1242c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1243c8b0795cSMark F. Adams   const PetscInt verbose      = pc_gamg->verbose;
1244c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1245c8b0795cSMark F. Adams   PetscInt       jj;
1246c5df96a5SBarry Smith   PetscMPIInt    rank,size;
1247c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
1248c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1249c8b0795cSMark F. Adams 
1250c8b0795cSMark F. Adams   PetscFunctionBegin;
12510cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12520cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
12530cbbd2e1SMark F. Adams #endif
1254c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
1255c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
1256c8b0795cSMark F. Adams 
12572e68589bSMark F. Adams   /* smooth P0 */
1258c8b0795cSMark F. Adams   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12592e68589bSMark F. Adams     Mat       tMat;
12602e68589bSMark F. Adams     Vec       diag;
12612e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12630cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12642e68589bSMark F. Adams #endif
12652e68589bSMark F. Adams     if (jj == 0) {
12662e68589bSMark F. Adams       KSP eksp;
12672e68589bSMark F. Adams       Vec bb, xx;
12682e68589bSMark F. Adams       PC  pc;
12692e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr);
12702e68589bSMark F. Adams       ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr);
12712e68589bSMark F. Adams       {
12722e68589bSMark F. Adams         PetscRandom rctx;
12732e68589bSMark F. Adams         ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
12742e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12752e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12762e68589bSMark F. Adams         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
12772e68589bSMark F. Adams       }
1278e696ed0bSMark F. Adams 
1279e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1280e696ed0bSMark F. Adams       {
1281e696ed0bSMark F. Adams         PetscInt    Istart,Iend,ncols,jj,Ii;
1282e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1283e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1284e696ed0bSMark F. Adams         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1285e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1286e696ed0bSMark F. Adams           if (ncols <= 1) {
1287e696ed0bSMark F. Adams             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1288e696ed0bSMark F. Adams           }
1289e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1290e696ed0bSMark F. Adams         }
1291e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1292e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1293e696ed0bSMark F. Adams       }
1294e696ed0bSMark F. Adams 
12952e68589bSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
1296806fa848SBarry Smith       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12972e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1298c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1299c2436cacSMark F. Adams       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1300c2436cacSMark F. Adams 
1301c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1302806fa848SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13032e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
13042e68589bSMark F. Adams 
1305c2436cacSMark F. Adams       ierr = KSPGetPC(eksp, &pc);CHKERRQ(ierr);
1306c2436cacSMark F. Adams       ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1307c2436cacSMark F. Adams 
13082e68589bSMark F. Adams       /* solve - keep stuff out of logging */
13092e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
13102e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
13112e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
13122e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
13132e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
13142e68589bSMark F. Adams 
13152e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
13169f3f12c8SMark F. Adams       if (verbose > 0) {
1317c8b0795cSMark F. Adams         PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
13182e68589bSMark F. Adams                     __FUNCT__,emax,emin,PCJACOBI);
13192e68589bSMark F. Adams       }
13202e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
13212e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
13222e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
13232e68589bSMark F. Adams 
13242e68589bSMark F. Adams       if (pc_gamg->emax_id == -1) {
13253bf036e2SBarry Smith         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
13262e68589bSMark F. Adams         assert(pc_gamg->emax_id != -1);
13272e68589bSMark F. Adams       }
13282e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
13292e68589bSMark F. Adams     }
13302e68589bSMark F. Adams 
13312e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
13322e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
13332e68589bSMark F. Adams     ierr  = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr);
13342e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
13352e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
13362e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
13372e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1338b4ec6429SMark F. Adams     alpha = -1.4/emax;
13392e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
13402e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
13412e68589bSMark F. Adams     Prol  = tMat;
13420cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13430cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
13442e68589bSMark F. Adams #endif
13452e68589bSMark F. Adams   }
13460cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13470cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
13480cbbd2e1SMark F. Adams #endif
1349c8b0795cSMark F. Adams   *a_P = Prol;
1350c8b0795cSMark F. Adams   PetscFunctionReturn(0);
13512e68589bSMark F. Adams }
13522e68589bSMark F. Adams 
1353c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1354c8b0795cSMark F. Adams /*
1355a2f3521dSMark F. Adams    PCGAMGKKTProl_AGG
1356a2f3521dSMark F. Adams 
1357a2f3521dSMark F. Adams   Input Parameter:
1358a2f3521dSMark F. Adams    . pc - this
1359a2f3521dSMark F. Adams    . Prol11 - matrix on this fine level
1360a2f3521dSMark F. Adams    . A21 - matrix on this fine level
1361a2f3521dSMark F. Adams  In/Output Parameter:
1362a2f3521dSMark F. Adams    . a_P22 - prolongation operator to the next level
1363a2f3521dSMark F. Adams */
1364a2f3521dSMark F. Adams #undef __FUNCT__
1365a2f3521dSMark F. Adams #define __FUNCT__ "PCGAMGKKTProl_AGG"
1366*2fa5cd67SKarl Rupp PetscErrorCode PCGAMGKKTProl_AGG(PC pc,const Mat Prol11,const Mat A21,Mat *a_P22)
1367a2f3521dSMark F. Adams {
1368a2f3521dSMark F. Adams   PetscErrorCode ierr;
1369a2f3521dSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1370a2f3521dSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1371a2f3521dSMark F. Adams   const PetscInt verbose  = pc_gamg->verbose;
1372a2f3521dSMark F. Adams   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1373c5df96a5SBarry Smith   PetscMPIInt      rank,size;
1374a2f3521dSMark F. Adams   Mat              Prol22,Tmat,Gmat;
1375a2f3521dSMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1376a2f3521dSMark F. Adams   PetscCoarsenData *agg_lists;
1377a2f3521dSMark F. Adams 
1378a2f3521dSMark F. Adams   PetscFunctionBegin;
1379a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1380a2f3521dSMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1381a2f3521dSMark F. Adams #endif
1382c5df96a5SBarry Smith   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
1383c5df96a5SBarry Smith   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
1384a2f3521dSMark F. Adams 
1385a2f3521dSMark F. Adams   /* form C graph */
1386a2f3521dSMark F. Adams   ierr = MatMatMult(A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);CHKERRQ(ierr);
1387a2f3521dSMark F. Adams   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat);CHKERRQ(ierr);
1388a2f3521dSMark F. Adams   ierr = MatDestroy(&Tmat);CHKERRQ(ierr);
1389a2f3521dSMark F. Adams   ierr = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose);CHKERRQ(ierr);
1390a2f3521dSMark F. Adams 
1391a2f3521dSMark F. Adams   /* coarsen constraints */
1392a2f3521dSMark F. Adams   {
1393a2f3521dSMark F. Adams     MatCoarsen crs;
1394a2f3521dSMark F. Adams     ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
1395a2f3521dSMark F. Adams     ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
1396a2f3521dSMark F. Adams     ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
1397a2f3521dSMark F. Adams     ierr = MatCoarsenSetVerbose(crs, verbose);CHKERRQ(ierr);
1398a2f3521dSMark F. Adams     ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1399a2f3521dSMark F. Adams     ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1400a2f3521dSMark F. Adams     ierr = MatCoarsenGetData(crs, &agg_lists);CHKERRQ(ierr);
1401a2f3521dSMark F. Adams     ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1402a2f3521dSMark F. Adams   }
1403a2f3521dSMark F. Adams 
1404a2f3521dSMark F. Adams   /* form simple prolongation 'Prol22' */
1405a2f3521dSMark F. Adams   {
1406a2f3521dSMark F. Adams     PetscInt    ii,mm,clid,my0,nloc,nLocalSelected;
1407a2f3521dSMark F. Adams     PetscScalar val = 1.0;
1408a2f3521dSMark F. Adams     /* get 'nLocalSelected' */
1409a2f3521dSMark F. Adams     ierr = MatGetLocalSize(Gmat, &nloc, &ii);CHKERRQ(ierr);
1410a2f3521dSMark F. Adams     for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1411a2f3521dSMark F. Adams       PetscBool ise;
1412a2f3521dSMark F. Adams       /* filter out singletons 0 or 1? */
1413a2f3521dSMark F. Adams       ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1414a2f3521dSMark F. Adams       if (!ise) nLocalSelected++;
1415a2f3521dSMark F. Adams     }
1416a2f3521dSMark F. Adams 
1417a2f3521dSMark F. Adams     ierr = MatCreate(wcomm,&Prol22);CHKERRQ(ierr);
1418806fa848SBarry Smith     ierr = MatSetSizes(Prol22,nloc, nLocalSelected,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1419a2f3521dSMark F. Adams     ierr = MatSetType(Prol22, MATAIJ);CHKERRQ(ierr);
1420a2f3521dSMark F. Adams     ierr = MatSeqAIJSetPreallocation(Prol22,1,PETSC_NULL);CHKERRQ(ierr);
1421806fa848SBarry Smith     ierr = MatMPIAIJSetPreallocation(Prol22,1,PETSC_NULL,1,PETSC_NULL);CHKERRQ(ierr);
1422a2f3521dSMark F. Adams     /* ierr = MatCreateAIJ(wcomm, */
1423a2f3521dSMark F. Adams     /*                      nloc, nLocalSelected, */
1424a2f3521dSMark F. Adams     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
1425a2f3521dSMark F. Adams     /*                      1, PETSC_NULL, 1, PETSC_NULL, */
1426a2f3521dSMark F. Adams     /*                      &Prol22); */
1427a2f3521dSMark F. Adams 
1428a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange(Prol22, &my0, &ii);CHKERRQ(ierr);
1429a2f3521dSMark F. Adams     nloc = ii - my0;
1430a2f3521dSMark F. Adams 
1431a2f3521dSMark F. Adams     /* make aggregates */
1432a2f3521dSMark F. Adams     for (mm = clid = 0; mm < nloc; mm++) {
1433a2f3521dSMark F. Adams       ierr = PetscCDSizeAt(agg_lists, mm, &ii);CHKERRQ(ierr);
1434a2f3521dSMark F. Adams       if (ii > 0) {
1435a2f3521dSMark F. Adams         PetscInt   asz=ii,cgid=my0+clid,rids[1000];
1436a2f3521dSMark F. Adams         PetscCDPos pos;
1437c666adbfSMark F. Adams         if (asz>1000) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1438a2f3521dSMark F. Adams         ii   = 0;
1439a2f3521dSMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1440a2f3521dSMark F. Adams         while (pos) {
1441a2f3521dSMark F. Adams           PetscInt gid1;
1442a2f3521dSMark F. Adams           ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
1443a2f3521dSMark F. Adams           ierr = PetscCDGetNextPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1444a2f3521dSMark F. Adams 
1445a2f3521dSMark F. Adams           rids[ii++] = gid1;
1446a2f3521dSMark F. Adams         }
1447a2f3521dSMark F. Adams         assert(ii==asz);
1448a2f3521dSMark F. Adams         /* add diagonal block of P0 */
1449a2f3521dSMark F. Adams         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES);CHKERRQ(ierr);
1450a2f3521dSMark F. Adams 
1451a2f3521dSMark F. Adams         clid++;
1452a2f3521dSMark F. Adams       } /* coarse agg */
1453a2f3521dSMark F. Adams     } /* for all fine nodes */
1454a2f3521dSMark F. Adams     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1455a2f3521dSMark F. Adams     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1456a2f3521dSMark F. Adams   }
1457a2f3521dSMark F. Adams 
1458a2f3521dSMark F. Adams   /* clean up */
1459a2f3521dSMark F. Adams   ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
1460a2f3521dSMark F. Adams   ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
1461a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1462a2f3521dSMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1463a2f3521dSMark F. Adams #endif
1464a2f3521dSMark F. Adams   *a_P22 = Prol22;
1465a2f3521dSMark F. Adams   PetscFunctionReturn(0);
1466a2f3521dSMark F. Adams }
1467a2f3521dSMark F. Adams 
1468a2f3521dSMark F. Adams /* -------------------------------------------------------------------------- */
1469a2f3521dSMark F. Adams /*
1470c8b0795cSMark F. Adams    PCCreateGAMG_AGG
14712e68589bSMark F. Adams 
1472c8b0795cSMark F. Adams   Input Parameter:
1473c8b0795cSMark F. Adams    . pc -
1474c8b0795cSMark F. Adams */
1475c8b0795cSMark F. Adams #undef __FUNCT__
1476c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1477c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1478c8b0795cSMark F. Adams {
1479c8b0795cSMark F. Adams   PetscErrorCode ierr;
1480c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1481c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1482c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
14832e68589bSMark F. Adams 
1484c8b0795cSMark F. Adams   PetscFunctionBegin;
1485ef820a72SMark F. Adams   if (pc_gamg->subctx) {
1486ef820a72SMark F. Adams     /* call base class */
1487ef820a72SMark F. Adams     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
1488ef820a72SMark F. Adams   }
1489ef820a72SMark F. Adams 
1490c8b0795cSMark F. Adams   /* create sub context for SA */
1491c8b0795cSMark F. Adams   ierr            = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr);
1492c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1493c8b0795cSMark F. Adams 
1494c8b0795cSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1495c8b0795cSMark F. Adams   pc->ops->destroy        = PCDestroy_AGG;
1496c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1497c8b0795cSMark F. Adams 
1498c8b0795cSMark F. Adams   /* set internal function pointers */
1499c8b0795cSMark F. Adams   pc_gamg->graph       = PCGAMGgraph_AGG;
1500b43b03e9SMark F. Adams   pc_gamg->coarsen     = PCGAMGCoarsen_AGG;
15010cbbd2e1SMark F. Adams   pc_gamg->prolongator = PCGAMGProlongator_AGG;
15020cbbd2e1SMark F. Adams   pc_gamg->optprol     = PCGAMGOptprol_AGG;
1503a2f3521dSMark F. Adams   pc_gamg->formkktprol = PCGAMGKKTProl_AGG;
1504c8b0795cSMark F. Adams 
1505c8b0795cSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_AGG;
1506c8b0795cSMark F. Adams 
15073bf036e2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_AGG",PCSetCoordinates_AGG);CHKERRQ(ierr);
15082e68589bSMark F. Adams   PetscFunctionReturn(0);
15092e68589bSMark F. Adams }
1510