xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 2d7fac45b48036e48b28675d0eaecfbd69e5e488)
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;
1742e68589bSMark F. Adams     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
1752e68589bSMark F. Adams                            "smoothing steps for smoothed aggregation, usually 1 (0)",
1762e68589bSMark F. Adams                            "PCGAMGSetNSmooths",
177c8b0795cSMark F. Adams                            pc_gamg_agg->nsmooths,
178c8b0795cSMark F. Adams                            &pc_gamg_agg->nsmooths,
179c8b0795cSMark F. Adams                            &flag);
180c8b0795cSMark F. Adams     CHKERRQ(ierr);
181c8b0795cSMark F. Adams 
182c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
183c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
184c8b0795cSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_sym_graph",
185581a99e3SJed Brown                             "Set for asymmetric matrices",
186c8b0795cSMark F. Adams                             "PCGAMGSetSymGraph",
187c8b0795cSMark F. Adams                             pc_gamg_agg->sym_graph,
188c8b0795cSMark F. Adams                             &pc_gamg_agg->sym_graph,
1892e68589bSMark F. Adams                             &flag);
1902e68589bSMark F. Adams     CHKERRQ(ierr);
191ef4ad70eSMark F. Adams 
192ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
193ef4ad70eSMark F. Adams     pc_gamg_agg->square_graph = PETSC_TRUE;
194ef4ad70eSMark F. Adams     ierr = PetscOptionsBool("-pc_gamg_square_graph",
1950cbbd2e1SMark F. Adams                             "For faster coarsening and lower coarse grid complexity",
196ef4ad70eSMark F. Adams                             "PCGAMGSetSquareGraph",
197ef4ad70eSMark F. Adams                             pc_gamg_agg->square_graph,
198ef4ad70eSMark F. Adams                             &pc_gamg_agg->square_graph,
199ef4ad70eSMark F. Adams                             &flag);
200ef4ad70eSMark F. Adams     CHKERRQ(ierr);
2012e68589bSMark F. Adams   }
2022e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
2032e68589bSMark F. Adams 
2042e68589bSMark F. Adams   PetscFunctionReturn(0);
2052e68589bSMark F. Adams }
2062e68589bSMark F. Adams 
2072e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2082e68589bSMark F. Adams /*
2092e68589bSMark F. Adams    PCDestroy_AGG
2102e68589bSMark F. Adams 
2112e68589bSMark F. Adams   Input Parameter:
2122e68589bSMark F. Adams    . pc -
2132e68589bSMark F. Adams */
2142e68589bSMark F. Adams #undef __FUNCT__
2152e68589bSMark F. Adams #define __FUNCT__ "PCDestroy_AGG"
2162e68589bSMark F. Adams PetscErrorCode PCDestroy_AGG( PC pc )
2172e68589bSMark F. Adams {
2182e68589bSMark F. Adams   PetscErrorCode  ierr;
2192e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
2202e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
221c8b0795cSMark F. Adams   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
2222e68589bSMark F. Adams 
2232e68589bSMark F. Adams   PetscFunctionBegin;
224c8b0795cSMark F. Adams   if( pc_gamg_agg ) {
225c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg_agg);CHKERRQ(ierr);
226c8b0795cSMark F. Adams     pc_gamg_agg = 0;
2272e68589bSMark F. Adams   }
2282e68589bSMark F. Adams 
2292e68589bSMark F. Adams   /* call base class */
2302e68589bSMark F. Adams   ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr);
2312e68589bSMark F. Adams 
2322e68589bSMark F. Adams   PetscFunctionReturn(0);
2332e68589bSMark F. Adams }
2342e68589bSMark F. Adams 
2352e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2362e68589bSMark F. Adams /*
2372e68589bSMark F. Adams    PCSetCoordinates_AGG
2382e68589bSMark F. Adams 
2392e68589bSMark F. Adams    Input Parameter:
2402e68589bSMark F. Adams    .  pc - the preconditioner context
2412e68589bSMark F. Adams */
2422e68589bSMark F. Adams EXTERN_C_BEGIN
2432e68589bSMark F. Adams #undef __FUNCT__
2442e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2452e68589bSMark F. Adams PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscReal *coords )
2462e68589bSMark F. Adams {
2472e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
2482e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2492e68589bSMark F. Adams   PetscErrorCode ierr;
2502e68589bSMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
2512e68589bSMark F. Adams   Mat            Amat = pc->pmat;
2522e68589bSMark F. Adams 
2532e68589bSMark F. Adams   PetscFunctionBegin;
2542e68589bSMark F. Adams   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
2552e68589bSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
2562e68589bSMark F. Adams   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
2572e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
2582e68589bSMark F. Adams   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
2592e68589bSMark F. Adams 
2602e68589bSMark F. Adams   /* SA: null space vectors */
261c8b0795cSMark F. Adams   if( coords && bs==1 ) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
262c8b0795cSMark F. Adams   else if( coords ) pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* elasticity */
263c8b0795cSMark F. Adams   else pc_gamg->data_cell_cols = bs; /* no data, force SA with constant null space vectors */
264c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = bs;
2652e68589bSMark F. Adams 
266c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2672e68589bSMark F. Adams 
2682e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2692e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2702e68589bSMark F. Adams     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
271c8b0795cSMark F. Adams     ierr = PetscMalloc(arrsz*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
2722e68589bSMark F. Adams   }
2732e68589bSMark F. Adams   /* copy data in - column oriented */
2742e68589bSMark F. Adams   for(kk=0;kk<nloc;kk++){
2752e68589bSMark F. Adams     const PetscInt M = Iend - my0;
2762e68589bSMark F. Adams     PetscReal *data = &pc_gamg->data[kk*bs];
277c8b0795cSMark F. Adams     if( pc_gamg->data_cell_cols==1 ) *data = 1.0;
2782e68589bSMark F. Adams     else {
2792e68589bSMark F. Adams       for(ii=0;ii<bs;ii++)
2802e68589bSMark F. Adams         for(jj=0;jj<bs;jj++)
2812e68589bSMark F. Adams           if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
2822e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2832e68589bSMark F. Adams       if( coords ) {
2842e68589bSMark F. Adams         if( ndm == 2 ){ /* rotational modes */
2852e68589bSMark F. Adams           data += 2*M;
2862e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2872e68589bSMark F. Adams           data[1] =  coords[2*kk];
2882e68589bSMark F. Adams         }
2892e68589bSMark F. Adams         else {
2902e68589bSMark F. Adams           data += 3*M;
2912e68589bSMark F. Adams           data[0] = 0.0;               data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2922e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  coords[3*kk];
2932e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2942e68589bSMark F. Adams         }
2952e68589bSMark F. Adams       }
2962e68589bSMark F. Adams     }
2972e68589bSMark F. Adams   }
2982e68589bSMark F. Adams 
2992e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
3002e68589bSMark F. Adams 
3012e68589bSMark F. Adams   PetscFunctionReturn(0);
3022e68589bSMark F. Adams }
3032e68589bSMark F. Adams EXTERN_C_END
3042e68589bSMark F. Adams 
305b43b03e9SMark F. Adams typedef PetscInt NState;
306b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
307b43b03e9SMark F. Adams static const NState DELETED=-1;
308b43b03e9SMark F. Adams static const NState REMOVED=-3;
309b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
310b43b03e9SMark F. Adams 
311c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
312c8b0795cSMark F. Adams /*
313b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
314b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
315c8b0795cSMark F. Adams 
316c8b0795cSMark F. Adams    Input Parameter:
317c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
318c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
319c8b0795cSMark F. Adams    Input/Output Parameter:
3200cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids )
321c8b0795cSMark F. Adams */
322c8b0795cSMark F. Adams #undef __FUNCT__
323c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3240cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */
3250cbbd2e1SMark F. Adams                                   const Mat Gmat_1, /* base graph */
3260cbbd2e1SMark F. Adams                                   /* const IS selected_2, [nselected local] selected vertices */
3270cbbd2e1SMark F. Adams                                   PetscCoarsenData *aggs_2 /* [nselected local] global ID of aggregate */
328c8b0795cSMark F. Adams                                   )
329c8b0795cSMark F. Adams {
330c8b0795cSMark F. Adams   PetscErrorCode ierr;
331c8b0795cSMark F. Adams   PetscBool      isMPI;
332c8b0795cSMark F. Adams   Mat_SeqAIJ    *matA_1, *matB_1=0, *matA_2, *matB_2=0;
333c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat_2)->comm;
3340cbbd2e1SMark F. Adams   PetscMPIInt    mype,npe;
3350cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
336c8b0795cSMark F. Adams   Mat_MPIAIJ    *mpimat_2 = 0, *mpimat_1=0;
337c8b0795cSMark F. Adams   const PetscInt nloc = Gmat_2->rmap->n;
3380cbbd2e1SMark F. Adams   PetscScalar   *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3390cbbd2e1SMark F. Adams   PetscInt      *lid_cprowID_1;
340c8b0795cSMark F. Adams   NState        *lid_state;
3410cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
342c8b0795cSMark F. Adams 
343c8b0795cSMark F. Adams   PetscFunctionBegin;
344c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype );   CHKERRQ(ierr);
3450cbbd2e1SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );   CHKERRQ(ierr);
346c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);  CHKERRQ(ierr);
347c8b0795cSMark F. Adams 
3480cbbd2e1SMark F. Adams   if( PETSC_FALSE ) {
349c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
350c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
351c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
352c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
353c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr);
354c8b0795cSMark F. Adams     ierr = PetscViewerDestroy( &viewer );
355c8b0795cSMark F. Adams   }
356c8b0795cSMark F. Adams 
357c8b0795cSMark F. Adams   /* get submatrices */
358c8b0795cSMark F. Adams   ierr = PetscTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr);
359c8b0795cSMark F. Adams   if(isMPI) {
360c8b0795cSMark F. Adams     /* grab matrix objects */
361c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
362c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
363c8b0795cSMark F. Adams     matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data;
364c8b0795cSMark F. Adams     matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data;
365c8b0795cSMark F. Adams     matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data;
366c8b0795cSMark F. Adams     matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data;
367c8b0795cSMark F. Adams 
368c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
369c8b0795cSMark F. Adams     matB_1->compressedrow.check = PETSC_TRUE;
370c8b0795cSMark F. Adams     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);
371c8b0795cSMark F. Adams     CHKERRQ(ierr);
372c8b0795cSMark F. Adams 
373c8b0795cSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr);
3740cbbd2e1SMark F. Adams     for( lid = 0 ; lid < nloc ; lid++ ) lid_cprowID_1[lid] = -1;
375c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
376c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
377c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
378c8b0795cSMark F. Adams     }
379c8b0795cSMark F. Adams   }
380c8b0795cSMark F. Adams   else {
381c8b0795cSMark F. Adams     matA_1 = (Mat_SeqAIJ*)Gmat_1->data;
382c8b0795cSMark F. Adams     matA_2 = (Mat_SeqAIJ*)Gmat_2->data;
3830cbbd2e1SMark F. Adams     lid_cprowID_1 = PETSC_NULL;
384c8b0795cSMark F. Adams   }
385c8b0795cSMark F. Adams   assert( matA_1 && !matA_1->compressedrow.use );
386c8b0795cSMark F. Adams   assert( matB_1==0 || matB_1->compressedrow.use );
387c8b0795cSMark F. Adams   assert( matA_2 && !matA_2->compressedrow.use );
388c8b0795cSMark F. Adams   assert( matB_2==0 || matB_2->compressedrow.use );
389c8b0795cSMark F. Adams 
390c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
391c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr);
3920cbbd2e1SMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscScalar), &lid_parent_gid ); CHKERRQ(ierr);
393c8b0795cSMark F. Adams   for( lid = 0 ; lid < nloc ; lid++ ) {
3940cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
395c8b0795cSMark F. Adams     lid_state[lid] = DELETED;
396c8b0795cSMark F. Adams   }
3970cbbd2e1SMark F. Adams 
3980cbbd2e1SMark F. Adams   /* set lid_state */
3990cbbd2e1SMark F. Adams   for( lid = 0 ; lid < nloc ; lid++ ) {
40041b27cdeSMark F. Adams     PetscCDPos pos;
401e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
402e78576d6SMark F. Adams     if( pos ) {
403e78576d6SMark F. Adams       PetscInt gid1;
404e78576d6SMark F. Adams       ierr = LLNGetID( pos, &gid1 ); CHKERRQ(ierr); assert(gid1==lid+my0);
4050cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
406b43b03e9SMark F. Adams     }
407b43b03e9SMark F. Adams   }
4080cbbd2e1SMark F. Adams 
4090cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
410c8b0795cSMark F. Adams   for(lid=kk=0;lid<nloc;lid++){
411c8b0795cSMark F. Adams     NState state = lid_state[lid];
412c8b0795cSMark F. Adams     if( IS_SELECTED(state) ){
41341b27cdeSMark F. Adams       PetscCDPos pos;
414e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
415e78576d6SMark F. Adams       while(pos){
416e78576d6SMark F. Adams         PetscInt gid1;
417e78576d6SMark F. Adams         ierr = LLNGetID( pos, &gid1 ); CHKERRQ(ierr);
418e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr);
419e78576d6SMark F. Adams 
4200cbbd2e1SMark F. Adams         if( gid1 >= my0 && gid1 < Iend ){
4210cbbd2e1SMark F. Adams           lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
422c8b0795cSMark F. Adams         }
423c8b0795cSMark F. Adams       }
4240cbbd2e1SMark F. Adams     }
4250cbbd2e1SMark F. Adams   }
4260cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
427c8b0795cSMark F. Adams   if (isMPI) {
428c8b0795cSMark F. Adams     Vec          tempVec;
429c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
430c8b0795cSMark F. Adams     ierr = MatGetVecs( Gmat_1, &tempVec, 0 );         CHKERRQ(ierr);
431c8b0795cSMark F. Adams     for(kk=0,j=my0;kk<nloc;kk++,j++){
432c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
433c8b0795cSMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
434c8b0795cSMark F. Adams     }
435c8b0795cSMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
436c8b0795cSMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
437c8b0795cSMark F. Adams     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
438c8b0795cSMark F. Adams     CHKERRQ(ierr);
439c8b0795cSMark F. Adams     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
440c8b0795cSMark F. Adams     CHKERRQ(ierr);
441c8b0795cSMark F. Adams     ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
442c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
443c8b0795cSMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
444c8b0795cSMark F. Adams     CHKERRQ(ierr);
445c8b0795cSMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
446c8b0795cSMark F. Adams     CHKERRQ(ierr);
447c8b0795cSMark F. Adams     ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
4480cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4490cbbd2e1SMark F. Adams     for(kk=0,j=my0;kk<nloc;kk++,j++){
4500cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4510cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
4520cbbd2e1SMark F. Adams     }
4530cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
4540cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
4550cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghost_par_orig2 ); CHKERRQ(ierr);
4560cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);
4570cbbd2e1SMark F. Adams     CHKERRQ(ierr);
4580cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);
4590cbbd2e1SMark F. Adams     CHKERRQ(ierr);
4600cbbd2e1SMark F. Adams     ierr = VecGetArray( ghost_par_orig2, &cpcol_2_par_orig ); CHKERRQ(ierr);
4610cbbd2e1SMark F. Adams 
462c8b0795cSMark F. Adams     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
463c8b0795cSMark F. Adams   } /* ismpi */
464c8b0795cSMark F. Adams 
465c8b0795cSMark F. Adams   /* doit */
466c8b0795cSMark F. Adams   for(lid=0;lid<nloc;lid++){
467c8b0795cSMark F. Adams     NState state = lid_state[lid];
4680cbbd2e1SMark F. Adams     if( IS_SELECTED(state) ) {
4690cbbd2e1SMark F. Adams       /* steal locals */
470c8b0795cSMark F. Adams       ii = matA_1->i; n = ii[lid+1] - ii[lid];
471c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
472c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4730cbbd2e1SMark F. Adams         PetscInt lidj = idx[j], sgid;
474c8b0795cSMark F. Adams         NState statej = lid_state[lidj];
4750cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4760cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4770cbbd2e1SMark F. Adams           if( sgid >= my0 && sgid < Iend ){       /* I'm stealing this local from a local sgid */
4780cbbd2e1SMark F. Adams             PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0;
47941b27cdeSMark F. Adams             PetscCDPos pos,last=PETSC_NULL;
480c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
481e78576d6SMark F. Adams             /* for( pos=PetscCDGetHeadPos(aggs_2,slid) ; pos ; pos=PetscCDGetNextPos(aggs_2,slid,pos)){ */
482e78576d6SMark F. Adams             /*   PetscInt gid = LLNGetID(pos); */
483e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos); CHKERRQ(ierr);
484e78576d6SMark F. Adams             while(pos){
485e78576d6SMark F. Adams               PetscInt gid;
486e78576d6SMark F. Adams               ierr = LLNGetID( 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;
493c8b0795cSMark F. Adams               }
4940cbbd2e1SMark F. Adams               else last = pos;
495e78576d6SMark F. Adams 
496e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos); CHKERRQ(ierr);
497c8b0795cSMark F. Adams             }
498c8b0795cSMark F. Adams             if(hav!=1){
499c8b0795cSMark F. Adams               if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
500c8b0795cSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav);
501c8b0795cSMark F. Adams             }
502c8b0795cSMark F. Adams           }
5030cbbd2e1SMark F. Adams           else{            /* I'm stealing this local, owned by a ghost */
504c8b0795cSMark F. Adams             assert(sgid==-1);
50541b27cdeSMark F. Adams             ierr = PetscCDAppendID( aggs_2, lid, lidj+my0 );      CHKERRQ(ierr);
506c8b0795cSMark F. Adams           }
507c8b0795cSMark F. Adams         }
5080cbbd2e1SMark F. Adams       } /* local neighbors */
509c8b0795cSMark F. Adams     }
510c8b0795cSMark F. Adams     else if( state == DELETED && lid_cprowID_1 ) {
5110cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
512c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
513c8b0795cSMark F. Adams       if( (ix=lid_cprowID_1[lid]) != -1 ){
514c8b0795cSMark F. Adams         ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
515c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
516c8b0795cSMark F. Adams         for( j=0 ; j<n ; j++ ) {
517c8b0795cSMark F. Adams           PetscInt cpid = idx[j];
518c8b0795cSMark F. Adams           NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
519c8b0795cSMark F. Adams           if( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */
5200cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
5210cbbd2e1SMark F. Adams             if( sgidold>=my0 && sgidold<Iend ) { /* this was mine */
5220cbbd2e1SMark F. Adams               PetscInt hav=0,oldslidj=sgidold-my0;
52341b27cdeSMark F. Adams               PetscCDPos pos,last=PETSC_NULL;
5240cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
525e78576d6SMark F. Adams               /* for( pos=PetscCDGetHeadPos(aggs_2,oldslidj) ; pos ; pos=PetscCDGetNextPos(aggs_2,oldslidj,pos)){ */
526e78576d6SMark F. Adams               /*   PetscInt gid = LLNGetID(pos); */
527e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr);
528e78576d6SMark F. Adams               while( pos ) {
529e78576d6SMark F. Adams                 PetscInt gid;
530e78576d6SMark F. Adams                 ierr = LLNGetID( pos, &gid ); CHKERRQ(ierr);
5310cbbd2e1SMark F. Adams                 if( lid+my0 == gid ) {
5320cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
5330cbbd2e1SMark F. Adams                   assert(last);
53441b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode( aggs_2, oldslidj, last ); CHKERRQ(ierr);
5350cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
5360cbbd2e1SMark F. Adams                   hav = 1;
537c8b0795cSMark F. Adams                   break;
538c8b0795cSMark F. Adams                 }
5390cbbd2e1SMark F. Adams                 else last = pos;
540e78576d6SMark F. Adams 
541e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr);
542c8b0795cSMark F. Adams               }
543c8b0795cSMark F. Adams               if(hav!=1){
544c8b0795cSMark F. Adams                 if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
545c8b0795cSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav);
546c8b0795cSMark F. Adams               }
547c8b0795cSMark F. Adams             }
5480cbbd2e1SMark F. Adams             else {
5490cbbd2e1SMark F. Adams               /* ghosts remove this later */
5500cbbd2e1SMark F. Adams             }
551c8b0795cSMark F. Adams           }
552c8b0795cSMark F. Adams         }
553c8b0795cSMark F. Adams       }
554c8b0795cSMark F. Adams     } /* selected/deleted */
555c8b0795cSMark F. Adams   } /* node loop */
556c8b0795cSMark F. Adams 
557c8b0795cSMark F. Adams   if( isMPI ) {
5580cbbd2e1SMark F. Adams     PetscScalar *cpcol_2_parent,*cpcol_2_gid;
5590cbbd2e1SMark F. Adams     Vec          tempVec,ghostgids2,ghostparents2;
5600cbbd2e1SMark F. Adams     PetscInt     cpid,nghost_2;
5610cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
562c8b0795cSMark F. Adams 
5630cbbd2e1SMark F. Adams     ierr = VecGetSize( mpimat_2->lvec, &nghost_2 );   CHKERRQ(ierr);
564c8b0795cSMark F. Adams     ierr = MatGetVecs( Gmat_2, &tempVec, 0 );         CHKERRQ(ierr);
5650cbbd2e1SMark F. Adams 
5660cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
567c8b0795cSMark F. Adams     for(kk=0,j=my0;kk<nloc;kk++,j++){
5680cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES );  CHKERRQ(ierr);
569c8b0795cSMark F. Adams     }
570c8b0795cSMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
571c8b0795cSMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
5720cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghostparents2 ); CHKERRQ(ierr);
5730cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);
574c8b0795cSMark F. Adams     CHKERRQ(ierr);
5750cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);
576c8b0795cSMark F. Adams     CHKERRQ(ierr);
5770cbbd2e1SMark F. Adams     ierr = VecGetArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr);
5780cbbd2e1SMark F. Adams 
5790cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5800cbbd2e1SMark F. Adams     for(kk=0,j=my0;kk<nloc;kk++,j++){
5810cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5820cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
5830cbbd2e1SMark F. Adams     }
5840cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
5850cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
5860cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghostgids2 ); CHKERRQ(ierr);
5870cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);
5880cbbd2e1SMark F. Adams     CHKERRQ(ierr);
5890cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);
5900cbbd2e1SMark F. Adams     CHKERRQ(ierr);
5910cbbd2e1SMark F. Adams     ierr = VecGetArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr);
5920cbbd2e1SMark F. Adams 
593c8b0795cSMark F. Adams     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
594c8b0795cSMark F. Adams 
5950cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
5960cbbd2e1SMark F. Adams     ierr = GAMGTableCreate( 2*nghost_2, &gid_cpid ); CHKERRQ(ierr);
5970cbbd2e1SMark F. Adams     for( cpid = 0 ; cpid < nghost_2 ; cpid++ ) {
5980cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5990cbbd2e1SMark F. Adams       if( state==DELETED ) {
6000cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6010cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
6020cbbd2e1SMark F. Adams         if( sgid_old == -1 && sgid_new != -1 ) {
6030cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6040cbbd2e1SMark F. Adams           ierr = GAMGTableAdd( &gid_cpid, gid, cpid ); CHKERRQ(ierr);
6050cbbd2e1SMark F. Adams         }
6060cbbd2e1SMark F. Adams       }
6070cbbd2e1SMark F. Adams     }
608c8b0795cSMark F. Adams 
6090cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
610c8b0795cSMark F. Adams     for(lid=0;lid<nloc;lid++){
611c8b0795cSMark F. Adams       NState state = lid_state[lid];
612c8b0795cSMark F. Adams       if( IS_SELECTED(state) ){
61341b27cdeSMark F. Adams         PetscCDPos pos,last=PETSC_NULL;
614c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
615e78576d6SMark F. Adams         /* for( pos=PetscCDGetHeadPos(aggs_2,lid) ; pos ; pos=PetscCDGetNextPos(aggs_2,lid,pos)){ */
616e78576d6SMark F. Adams         /*   PetscInt gid = LLNGetID(pos); */
617e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
618e78576d6SMark F. Adams         while(pos){
619e78576d6SMark F. Adams           PetscInt gid;
620e78576d6SMark F. Adams           ierr = LLNGetID( pos, &gid ); CHKERRQ(ierr);
621e78576d6SMark F. Adams 
6220cbbd2e1SMark F. Adams           if( gid < my0 || gid >= Iend ) {
6230cbbd2e1SMark F. Adams             ierr = GAMGTableFind( &gid_cpid, gid, &cpid ); CHKERRQ(ierr);
6240cbbd2e1SMark F. Adams             if( cpid != -1 ) {
6250cbbd2e1SMark F. Adams               /* a moved ghost - */
6260cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
62741b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode( aggs_2, lid, last ); CHKERRQ(ierr);
6280cbbd2e1SMark F. Adams             }
6290cbbd2e1SMark F. Adams             else last = pos;
6300cbbd2e1SMark F. Adams           }
6310cbbd2e1SMark F. Adams           else last = pos;
632e78576d6SMark F. Adams 
633e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr);
634c8b0795cSMark F. Adams         } /* loop over list of deleted */
635c8b0795cSMark F. Adams       } /* selected */
636c8b0795cSMark F. Adams     }
6370cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy( &gid_cpid ); CHKERRQ(ierr);
638c8b0795cSMark F. Adams 
6390cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
6400cbbd2e1SMark F. Adams     for( cpid = 0 ; cpid < nghost_2 ; cpid++ ){
6410cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6420cbbd2e1SMark F. Adams       if( sgid_new >= my0 && sgid_new < Iend ) { /* this is mine */
6430cbbd2e1SMark F. Adams         PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6440cbbd2e1SMark F. Adams         PetscInt slid_new=sgid_new-my0,hav=0;
64541b27cdeSMark F. Adams         PetscCDPos pos;
6460cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
647e78576d6SMark F. Adams         /* for( pos=PetscCDGetHeadPos(aggs_2,slid_new) ; pos ; pos=PetscCDGetNextPos(aggs_2,slid_new,pos) ) { */
648e78576d6SMark F. Adams         /*   PetscInt gidj = LLNGetID(pos); */
649e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos); CHKERRQ(ierr);
650e78576d6SMark F. Adams         while(pos){
651e78576d6SMark F. Adams           PetscInt gidj;
652e78576d6SMark F. Adams           ierr = LLNGetID( pos, &gidj ); CHKERRQ(ierr);
653e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos); CHKERRQ(ierr);
654e78576d6SMark F. Adams 
6550cbbd2e1SMark F. Adams           if( gidj == gid ) { hav = 1; break; }
656c8b0795cSMark F. Adams         }
657c8b0795cSMark F. Adams         if( hav != 1 ){
6580cbbd2e1SMark F. Adams           /* id_llist_2[flidj] = id_llist_2[slid_new]; id_llist_2[slid_new] = flidj; /\* insert 'flidj' into head of llist *\/ */
65941b27cdeSMark F. Adams           ierr = PetscCDAppendID( aggs_2, slid_new, gid );      CHKERRQ(ierr);
660c8b0795cSMark F. Adams         }
661c8b0795cSMark F. Adams       }
662c8b0795cSMark F. Adams     }
663c8b0795cSMark F. Adams 
6640cbbd2e1SMark F. Adams     ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
6650cbbd2e1SMark F. Adams     ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
6660cbbd2e1SMark F. Adams     ierr = VecRestoreArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr);
6670cbbd2e1SMark F. Adams     ierr = VecRestoreArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr);
668c8b0795cSMark F. Adams     ierr = PetscFree( lid_cprowID_1 );  CHKERRQ(ierr);
6690cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghostgids2 ); CHKERRQ(ierr);
6700cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghostparents2 ); CHKERRQ(ierr);
6710cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghost_par_orig2 ); CHKERRQ(ierr);
672c8b0795cSMark F. Adams   }
673c8b0795cSMark F. Adams 
6740cbbd2e1SMark F. Adams   ierr = PetscFree( lid_parent_gid );  CHKERRQ(ierr);
675c8b0795cSMark F. Adams   ierr = PetscFree( lid_state );  CHKERRQ(ierr);
676c8b0795cSMark F. Adams 
677c8b0795cSMark F. Adams   PetscFunctionReturn(0);
678c8b0795cSMark F. Adams }
6792e68589bSMark F. Adams 
6802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6812e68589bSMark F. Adams /*
6822e68589bSMark F. Adams    PCSetData_AGG
6832e68589bSMark F. Adams 
6842e68589bSMark F. Adams   Input Parameter:
6852e68589bSMark F. Adams    . pc -
6862e68589bSMark F. Adams */
6872e68589bSMark F. Adams #undef __FUNCT__
6882e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
6892e68589bSMark F. Adams PetscErrorCode PCSetData_AGG( PC pc )
6902e68589bSMark F. Adams {
6912e68589bSMark F. Adams   PetscErrorCode  ierr;
6922e68589bSMark F. Adams   PetscFunctionBegin;
6932e68589bSMark F. Adams   ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr);
6942e68589bSMark F. Adams   PetscFunctionReturn(0);
6952e68589bSMark F. Adams }
6962e68589bSMark F. Adams 
6972e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6982e68589bSMark F. Adams /*
6992e68589bSMark F. Adams  formProl0
7002e68589bSMark F. Adams 
7012e68589bSMark F. Adams    Input Parameter:
7020cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
7032e68589bSMark F. Adams    . bs - block size
7040cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7050cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7062e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
7072e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7082e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7092e68589bSMark F. Adams   Output Parameter:
7102e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7112e68589bSMark F. Adams    . a_Prol - prolongation operator
7122e68589bSMark F. Adams */
7132e68589bSMark F. Adams #undef __FUNCT__
7142e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7150cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists,/* list from selected vertices of aggregate unselected vertices */
7160cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
7170cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
7180cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
7190cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7200cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7210cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7222e68589bSMark F. Adams                                 PetscReal **a_data_out,
7232e68589bSMark F. Adams                                 Mat a_Prol /* prolongation operator (output)*/
7242e68589bSMark F. Adams                                 )
7252e68589bSMark F. Adams {
7262e68589bSMark F. Adams   PetscErrorCode ierr;
7270cbbd2e1SMark F. Adams   PetscInt  Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7282e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
7292e68589bSMark F. Adams   PetscMPIInt    mype, npe;
7302e68589bSMark F. Adams   PetscReal      *out_data;
73141b27cdeSMark F. Adams   PetscCDPos         pos;
7320cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7330cbbd2e1SMark F. Adams 
7342adcac29SMark F. Adams /* #define OUT_AGGS */
7359057884aSMark F. Adams #ifdef OUT_AGGS
7360cbbd2e1SMark F. Adams   static PetscInt llev = 0; char fname[32]; FILE *file; PetscInt pM;
7379057884aSMark F. Adams #endif
7382e68589bSMark F. Adams 
7392e68589bSMark F. Adams   PetscFunctionBegin;
7402e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
7412e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
7422e68589bSMark F. Adams   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );    CHKERRQ(ierr);
7430cbbd2e1SMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
7440cbbd2e1SMark F. Adams   Iend /= bs;
7450cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7462e68589bSMark F. Adams 
7470cbbd2e1SMark F. Adams   ierr = GAMGTableCreate( 2*nghosts, &fgid_flid ); CHKERRQ(ierr);
7480cbbd2e1SMark F. Adams   for(kk=0;kk<nghosts;kk++) {
7490cbbd2e1SMark F. Adams     ierr = GAMGTableAdd( &fgid_flid, flid_fgid[nloc+kk], nloc+kk ); CHKERRQ(ierr);
7502e68589bSMark F. Adams   }
7512e68589bSMark F. Adams 
7520cbbd2e1SMark F. Adams #ifdef OUT_AGGS
7530cbbd2e1SMark F. Adams   sprintf(fname,"aggs_%d_%d.m",llev++,mype);
7540cbbd2e1SMark F. Adams   if(llev==1) {
7550cbbd2e1SMark F. Adams     file = fopen(fname,"w");
7560cbbd2e1SMark F. Adams   }
7570cbbd2e1SMark F. Adams   MatGetSize( a_Prol, &pM, &jj );
7580cbbd2e1SMark F. Adams #endif
7590cbbd2e1SMark F. Adams 
7600cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7610cbbd2e1SMark F. Adams   for(nSelected=mm=0;mm<nloc;mm++) {
762e78576d6SMark F. Adams     PetscBool ise;
763e78576d6SMark F. Adams     ierr = PetscCDEmptyAt( agg_llists, mm, &ise ); CHKERRQ(ierr);
764e78576d6SMark F. Adams     if( !ise ) nSelected++;
7650cbbd2e1SMark F. Adams   }
7660cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn( a_Prol, &ii, &jj ); CHKERRQ(ierr);
7670cbbd2e1SMark F. Adams   assert((ii/nSAvec)==my0crs); assert(nSelected==(jj-ii)/nSAvec);
7680cbbd2e1SMark F. Adams 
7692e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7700cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7710cbbd2e1SMark F. Adams   ierr = PetscMalloc( out_data_stride*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr);
7720cbbd2e1SMark F. Adams   for(ii=0;ii<out_data_stride*nSAvec;ii++) {
7730cbbd2e1SMark F. Adams     out_data[ii]=1.e300;
7740cbbd2e1SMark F. Adams   }
7750cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7762e68589bSMark F. Adams 
7772e68589bSMark F. Adams   /* find points and set prolongation */
778c8b0795cSMark F. Adams   minsz = 100;
7792e68589bSMark F. Adams   ndone = 0;
7800cbbd2e1SMark F. Adams   for( mm = clid = 0 ; mm < nloc ; mm++ ){
781e78576d6SMark F. Adams     ierr = PetscCDSizeAt( agg_llists, mm, &jj ); CHKERRQ(ierr);
782e78576d6SMark F. Adams     if( jj > 0 ) {
7830cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7840cbbd2e1SMark F. Adams       PetscInt cids[100]; /* max bs */
7850cbbd2e1SMark F. Adams       PetscBLASInt asz=jj,M=asz*bs,N=nSAvec,INFO;
7862e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs;
7872e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7882e68589bSMark F. Adams       PetscInt       *fids;
78965d7b583SSatish Balay       PetscReal      *data;
7900cbbd2e1SMark F. Adams       /* count agg */
7910cbbd2e1SMark F. Adams       if( asz<minsz ) minsz = asz;
7920cbbd2e1SMark F. Adams 
7930cbbd2e1SMark F. Adams       /* get block */
7942e68589bSMark F. Adams       ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr);
7952e68589bSMark F. Adams       ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr);
7962e68589bSMark F. Adams       ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr);
7972e68589bSMark F. Adams       ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr);
7982e68589bSMark F. Adams       ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr);
7992e68589bSMark F. Adams 
8002e68589bSMark F. Adams       aggID = 0;
801e78576d6SMark F. Adams       /* for( pos=PetscCDGetHeadPos(agg_llists,lid) ;  */
802e78576d6SMark F. Adams       /*      pos ;  */
803e78576d6SMark F. Adams       /*      pos=PetscCDGetNextPos(agg_llists,lid,pos)) { */
804e78576d6SMark F. Adams       /*   PetscInt gid1 = LLNGetID(pos); */
805e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(agg_llists,lid,&pos); CHKERRQ(ierr);
806e78576d6SMark F. Adams       while(pos){
807e78576d6SMark F. Adams         PetscInt gid1;
808e78576d6SMark F. Adams         ierr = LLNGetID( pos, &gid1 ); CHKERRQ(ierr);
809e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos); CHKERRQ(ierr);
810e78576d6SMark F. Adams 
8110cbbd2e1SMark F. Adams         if( gid1 >= my0 && gid1 < Iend ) flid = gid1 - my0;
8120cbbd2e1SMark F. Adams         else {
8130cbbd2e1SMark F. Adams           ierr = GAMGTableFind( &fgid_flid, gid1, &flid ); CHKERRQ(ierr);
8140cbbd2e1SMark F. Adams           assert(flid>=0);
8150cbbd2e1SMark F. Adams         }
8162e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
81765d7b583SSatish Balay         data = &data_in[flid*bs];
8182e68589bSMark F. Adams         for( kk = ii = 0; ii < bs ; ii++ ) {
8192e68589bSMark F. Adams           for( jj = 0; jj < N ; jj++ ) {
8200cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8210cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8222e68589bSMark F. Adams           }
8232e68589bSMark F. Adams         }
8249057884aSMark F. Adams #ifdef OUT_AGGS
825b2a4f308SMark F. Adams         if(llev==1) {
8269057884aSMark F. Adams           char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8270cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
8280cbbd2e1SMark F. Adams           str[12] = col[(clid+17*mype)%6]; str[13] = sim[(clid+17*mype)%11];
8290cbbd2e1SMark F. Adams           MM = (PetscInt)(PetscSqrtScalar((PetscScalar)pM));
8300cbbd2e1SMark F. Adams           pj = gid1/MM; pi = gid1%MM;
831b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
832b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8339057884aSMark F. Adams         }
8349057884aSMark F. Adams #endif
8352e68589bSMark F. Adams         /* set fine IDs */
8362e68589bSMark F. Adams         for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8372e68589bSMark F. Adams 
8382e68589bSMark F. Adams         aggID++;
8390cbbd2e1SMark F. Adams       }
8402e68589bSMark F. Adams 
8412e68589bSMark F. Adams       /* pad with zeros */
8422e68589bSMark F. Adams       for( ii = asz*bs; ii < Mdata ; ii++ ) {
8432e68589bSMark F. Adams 	for( jj = 0; jj < N ; jj++, kk++ ) {
8442e68589bSMark F. Adams 	  qqc[jj*Mdata + ii] = .0;
8452e68589bSMark F. Adams 	}
8462e68589bSMark F. Adams       }
8472e68589bSMark F. Adams 
8482e68589bSMark F. Adams       ndone += aggID;
8492e68589bSMark F. Adams       /* QR */
8502e68589bSMark F. Adams       LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
8512e68589bSMark F. Adams       if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error");
8522e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8532e68589bSMark F. Adams       {
8542e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8552e68589bSMark F. Adams         for( jj = 0; jj < nSAvec ; jj++ ) {
8562e68589bSMark F. Adams           for( ii = 0; ii < nSAvec ; ii++ ) {
8570cbbd2e1SMark F. Adams             assert(data[jj*out_data_stride + ii] == 1.e300);
8580cbbd2e1SMark F. Adams             if( ii <= jj ) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8590cbbd2e1SMark F. Adams 	    else data[jj*out_data_stride + ii] = 0.;
8602e68589bSMark F. Adams           }
8612e68589bSMark F. Adams         }
8622e68589bSMark F. Adams       }
8632e68589bSMark F. Adams 
8642e68589bSMark F. Adams       /* get Q - row oriented */
8652e68589bSMark F. Adams       LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
8662e68589bSMark F. Adams       if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO);
8672e68589bSMark F. Adams 
8682e68589bSMark F. Adams       for( ii = 0 ; ii < M ; ii++ ){
8692e68589bSMark F. Adams         for( jj = 0 ; jj < N ; jj++ ) {
8702e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8712e68589bSMark F. Adams         }
8722e68589bSMark F. Adams       }
8732e68589bSMark F. Adams 
8742e68589bSMark F. Adams       /* add diagonal block of P0 */
875c8b0795cSMark F. Adams       for(kk=0;kk<N;kk++) {
876c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
877c8b0795cSMark F. Adams       }
8782e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr);
8792e68589bSMark F. Adams 
8802e68589bSMark F. Adams       ierr = PetscFree( qqc );  CHKERRQ(ierr);
8812e68589bSMark F. Adams       ierr = PetscFree( qqr );  CHKERRQ(ierr);
8822e68589bSMark F. Adams       ierr = PetscFree( TAU );  CHKERRQ(ierr);
8832e68589bSMark F. Adams       ierr = PetscFree( WORK );  CHKERRQ(ierr);
8842e68589bSMark F. Adams       ierr = PetscFree( fids );  CHKERRQ(ierr);
885b43b03e9SMark F. Adams       clid++;
8860cbbd2e1SMark F. Adams     } /* coarse agg */
8870cbbd2e1SMark F. Adams   } /* for all fine nodes */
8880cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8890cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8902e68589bSMark F. Adams 
891c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */
8922e68589bSMark F. Adams /* MatGetSize( a_Prol, &kk, &jj ); */
893c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */
894e78576d6SMark F. Adams /* PetscPrintf(wcomm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",mype,__FUNCT__,ii,kk/bs,ndone,jj); */
8952e68589bSMark F. Adams 
8969057884aSMark F. Adams #ifdef OUT_AGGS
897b2a4f308SMark F. Adams   if(llev==1) fclose(file);
8989057884aSMark F. Adams #endif
8990cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy( &fgid_flid ); CHKERRQ(ierr);
9002e68589bSMark F. Adams 
9012e68589bSMark F. Adams   PetscFunctionReturn(0);
9022e68589bSMark F. Adams }
9032e68589bSMark F. Adams 
9042e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
9052e68589bSMark F. Adams /*
906c8b0795cSMark F. Adams    PCGAMGgraph_AGG
9072e68589bSMark F. Adams 
9082e68589bSMark F. Adams   Input Parameter:
9092e68589bSMark F. Adams    . pc - this
9102e68589bSMark F. Adams    . Amat - matrix on this fine level
9112e68589bSMark F. Adams   Output Parameter:
912c8b0795cSMark F. Adams    . a_Gmat -
9132e68589bSMark F. Adams */
9142e68589bSMark F. Adams #undef __FUNCT__
915c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG"
916c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_AGG( PC pc,
9172e68589bSMark F. Adams                                 const Mat Amat,
918c8b0795cSMark F. Adams                                 Mat *a_Gmat
919c8b0795cSMark F. Adams                                 )
920c8b0795cSMark F. Adams {
921c8b0795cSMark F. Adams   PetscErrorCode ierr;
922c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
923c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
924c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
925c8b0795cSMark F. Adams   const PetscReal vfilter = pc_gamg->threshold;
926c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
927c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
928e0940f08SMark F. Adams   Mat            Gmat;
929c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
9300cbbd2e1SMark F. Adams   PetscBool  set,flg,symm;
931c8b0795cSMark F. Adams 
932c8b0795cSMark F. Adams   PetscFunctionBegin;
9330cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9340cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9350cbbd2e1SMark F. Adams #endif
936c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
937c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
938c8b0795cSMark F. Adams 
9390cbbd2e1SMark F. Adams   ierr = MatIsSymmetricKnown(Amat, &set, &flg);        CHKERRQ(ierr);
940263489e9SJed Brown   symm = (PetscBool)(pc_gamg_agg->sym_graph || !(set && flg));
9410cbbd2e1SMark F. Adams 
942*2d7fac45SMark F. Adams   ierr  = PCGAMGCreateGraph( Amat, &Gmat ); CHKERRQ( ierr );
943*2d7fac45SMark F. Adams   ierr  = PCGAMGFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr );
944c8b0795cSMark F. Adams 
945e0940f08SMark F. Adams   *a_Gmat = Gmat;
946c8b0795cSMark F. Adams 
9470cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9480cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9490cbbd2e1SMark F. Adams #endif
950c8b0795cSMark F. Adams   PetscFunctionReturn(0);
951c8b0795cSMark F. Adams }
952c8b0795cSMark F. Adams 
953c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
954c8b0795cSMark F. Adams /*
955b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
956c8b0795cSMark F. Adams 
957c8b0795cSMark F. Adams   Input Parameter:
958e0940f08SMark F. Adams    . a_pc - this
959e0940f08SMark F. Adams   Input/Output Parameter:
9600cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
961c8b0795cSMark F. Adams   Output Parameter:
9620cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
963c8b0795cSMark F. Adams */
964c8b0795cSMark F. Adams #undef __FUNCT__
965b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
966e0940f08SMark F. Adams PetscErrorCode PCGAMGCoarsen_AGG( PC a_pc,
967e0940f08SMark F. Adams                                   Mat *a_Gmat1,
9680cbbd2e1SMark F. Adams                                   PetscCoarsenData **agg_lists
969c8b0795cSMark F. Adams                                   )
970c8b0795cSMark F. Adams {
971c8b0795cSMark F. Adams   PetscErrorCode ierr;
972e0940f08SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
973c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
974c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9750cbbd2e1SMark F. Adams   Mat             mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */
9760cbbd2e1SMark F. Adams   IS              perm;
977c8b0795cSMark F. Adams   PetscInt        Ii,nloc,bs,n,m;
978c8b0795cSMark F. Adams   PetscInt *permute;
979c8b0795cSMark F. Adams   PetscBool *bIndexSet;
980b43b03e9SMark F. Adams   MatCoarsen crs;
981e0940f08SMark F. Adams   MPI_Comm        wcomm = ((PetscObject)Gmat1)->comm;
9820cbbd2e1SMark F. Adams   PetscMPIInt     mype,npe;
983c8b0795cSMark F. Adams 
984c8b0795cSMark F. Adams   PetscFunctionBegin;
9850cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9860cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9870cbbd2e1SMark F. Adams #endif
9880cbbd2e1SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
9890cbbd2e1SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
990e0940f08SMark F. Adams   ierr = MatGetLocalSize( Gmat1, &n, &m ); CHKERRQ(ierr);
991e0940f08SMark F. Adams   ierr = MatGetBlockSize( Gmat1, &bs ); CHKERRQ(ierr); assert(bs==1);
992c8b0795cSMark F. Adams   nloc = n/bs;
993c8b0795cSMark F. Adams 
994e0940f08SMark F. Adams   if( pc_gamg_agg->square_graph ) {
995e0940f08SMark F. Adams     ierr = MatTransposeMatMult( Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 );
996e0940f08SMark F. Adams     CHKERRQ(ierr);
997e0940f08SMark F. Adams   }
998e0940f08SMark F. Adams   else Gmat2 = Gmat1;
999c8b0795cSMark F. Adams 
1000c8b0795cSMark F. Adams   /* get MIS aggs */
1001c8b0795cSMark F. Adams   /* randomize */
1002c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
1003c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
1004c8b0795cSMark F. Adams   for ( Ii = 0; Ii < nloc ; Ii++ ){
1005c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
1006c8b0795cSMark F. Adams     permute[Ii] = Ii;
1007c8b0795cSMark F. Adams   }
1008c8b0795cSMark F. Adams   srand(1); /* make deterministic */
1009c8b0795cSMark F. Adams   for ( Ii = 0; Ii < nloc ; Ii++ ) {
1010c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
1011c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1012c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1013c8b0795cSMark F. Adams       permute[iSwapIndex] = permute[Ii];
1014c8b0795cSMark F. Adams       permute[Ii] = iTemp;
1015c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1016c8b0795cSMark F. Adams     }
1017c8b0795cSMark F. Adams   }
1018c8b0795cSMark F. Adams   ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
1019c8b0795cSMark F. Adams 
1020c8b0795cSMark F. Adams   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);
1021c8b0795cSMark F. Adams   CHKERRQ(ierr);
10220cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10230cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1024b43b03e9SMark F. Adams #endif
1025b43b03e9SMark F. Adams   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
10269057884aSMark F. Adams   /* ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); */
10279057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions( crs ); CHKERRQ(ierr);
1028b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
1029b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr);
1030b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
1031b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr);
1032b43b03e9SMark F. Adams   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
10330cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData( crs, agg_lists ); CHKERRQ(ierr); /* output */
1034b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
1035b43b03e9SMark F. Adams 
1036c8b0795cSMark F. Adams   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
1037c8b0795cSMark F. Adams   ierr = PetscFree( permute );  CHKERRQ(ierr);
10380cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10390cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1040b43b03e9SMark F. Adams #endif
1041c8b0795cSMark F. Adams   /* smooth aggs */
1042e0940f08SMark F. Adams   if( Gmat2 != Gmat1 ) {
10430cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10440cbbd2e1SMark F. Adams     ierr = smoothAggs( Gmat2, Gmat1, *agg_lists ); CHKERRQ(ierr);
1045c8b0795cSMark F. Adams     ierr = MatDestroy( &Gmat1 );  CHKERRQ(ierr);
1046e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
104741b27cdeSMark F. Adams     ierr = PetscCDGetMat( llist, &mat );  CHKERRQ(ierr);
10480cbbd2e1SMark F. Adams     if(mat) SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1049ef4ad70eSMark F. Adams   }
10500cbbd2e1SMark F. Adams   else {
10510cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10520cbbd2e1SMark F. Adams     /* see if we have a matrix that takes pecedence (returned from MatCoarsenAppply) */
105341b27cdeSMark F. Adams     ierr = PetscCDGetMat( llist, &mat );   CHKERRQ(ierr);
10540cbbd2e1SMark F. Adams     if( mat ) {
10550cbbd2e1SMark F. Adams       ierr = MatDestroy( &Gmat1 );  CHKERRQ(ierr);
10560cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10570cbbd2e1SMark F. Adams     }
10580cbbd2e1SMark F. Adams   }
10590cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
10600cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
10610cbbd2e1SMark F. Adams #endif
1062c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1063c8b0795cSMark F. Adams }
1064c8b0795cSMark F. Adams 
1065c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1066c8b0795cSMark F. Adams /*
10670cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1068c8b0795cSMark F. Adams 
1069c8b0795cSMark F. Adams  Input Parameter:
1070c8b0795cSMark F. Adams  . pc - this
1071c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1072c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10730cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1074c8b0795cSMark F. Adams  Output Parameter:
1075c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1076c8b0795cSMark F. Adams  */
1077c8b0795cSMark F. Adams #undef __FUNCT__
10780cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
10790cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_AGG( PC pc,
1080c8b0795cSMark F. Adams                                       const Mat Amat,
1081c8b0795cSMark F. Adams                                       const Mat Gmat,
10820cbbd2e1SMark F. Adams                                       PetscCoarsenData *agg_lists,
1083c8b0795cSMark F. Adams                                       Mat *a_P_out
10842e68589bSMark F. Adams                                       )
10852e68589bSMark F. Adams {
10862e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
10872e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
10882e68589bSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
1089c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
10902e68589bSMark F. Adams   PetscErrorCode ierr;
1091c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1092c8b0795cSMark F. Adams   Mat            Prol;
10932e68589bSMark F. Adams   PetscMPIInt    mype, npe;
10942e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
10950cbbd2e1SMark F. Adams   const PetscInt col_bs=data_cols;
1096c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1097c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
10982e68589bSMark F. Adams 
10992e68589bSMark F. Adams   PetscFunctionBegin;
11000cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11010cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11020cbbd2e1SMark F. Adams #endif
11032e68589bSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
11042e68589bSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
11052e68589bSMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
1106c8b0795cSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr );
1107c8b0795cSMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
11082e68589bSMark F. Adams 
11092e68589bSMark F. Adams   /* get 'nLocalSelected' */
11100cbbd2e1SMark F. Adams   for( ii=0, nLocalSelected = 0 ; ii < nloc ; ii++ ){
1111e78576d6SMark F. Adams     PetscBool ise;
11120cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1113e78576d6SMark F. Adams     ierr = PetscCDEmptyAt( agg_lists, ii, &ise ); CHKERRQ(ierr);
1114e78576d6SMark F. Adams     if( !ise ) nLocalSelected++;
11152e68589bSMark F. Adams   }
11162e68589bSMark F. Adams 
11172e68589bSMark F. Adams   /* create prolongator, create P matrix */
111869b1f4b7SBarry Smith   ierr = MatCreateAIJ( wcomm,
1119c8b0795cSMark F. Adams                        nloc*bs, nLocalSelected*col_bs,
11202e68589bSMark F. Adams                        PETSC_DETERMINE, PETSC_DETERMINE,
11212e68589bSMark F. Adams                        data_cols, PETSC_NULL, data_cols, PETSC_NULL,
11222e68589bSMark F. Adams                        &Prol );
11232e68589bSMark F. Adams   CHKERRQ(ierr);
11242e68589bSMark F. Adams 
11252e68589bSMark F. Adams   /* can get all points "removed" */
1126c8b0795cSMark F. Adams   ierr =  MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr);
1127c8b0795cSMark F. Adams   if( ii==0 ) {
11282e68589bSMark F. Adams     if( verbose ) {
1129c8b0795cSMark F. Adams       PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__);
11302e68589bSMark F. Adams     }
11312e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
11322e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
11332e68589bSMark F. Adams     PetscFunctionReturn(0);
11342e68589bSMark F. Adams   }
1135c8b0795cSMark F. Adams   if( verbose ) {
1136e78576d6SMark F. Adams     PetscPrintf(wcomm,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/col_bs);
1137c8b0795cSMark F. Adams   }
1138c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr);
11390cbbd2e1SMark F. Adams 
11400cbbd2e1SMark F. Adams   assert((kk-myCrs0)%col_bs==0);
1141c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
11420cbbd2e1SMark F. Adams   assert((kk/col_bs-myCrs0)==nLocalSelected);
11432e68589bSMark F. Adams 
11442e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11460cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11472e68589bSMark F. Adams #endif
1148c8b0795cSMark F. Adams   if (npe > 1) { /*  */
11492e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
11502e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr);
11512e68589bSMark F. Adams     for( jj = 0 ; jj < data_cols ; jj++ ){
1152c8b0795cSMark F. Adams       for( kk = 0 ; kk < bs ; kk++) {
11532e68589bSMark F. Adams         PetscInt ii,nnodes;
1154c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1155c8b0795cSMark F. Adams         for( ii = 0 ; ii < nloc ; ii++, tp += bs ){
11562e68589bSMark F. Adams           tmp_ldata[ii] = *tp;
11572e68589bSMark F. Adams         }
11580cbbd2e1SMark F. Adams         ierr = PCGAMGGetDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata );
11592e68589bSMark F. Adams         CHKERRQ(ierr);
11602e68589bSMark F. Adams         if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1161c8b0795cSMark F. Adams           ierr = PetscMalloc( nnodes*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr);
1162c8b0795cSMark F. Adams           nbnodes = bs*nnodes;
11632e68589bSMark F. Adams         }
1164c8b0795cSMark F. Adams         tp2 = data_w_ghost + jj*bs*nnodes + kk;
1165c8b0795cSMark F. Adams         for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs ){
11662e68589bSMark F. Adams           *tp2 = tmp_gdata[ii];
11672e68589bSMark F. Adams         }
11682e68589bSMark F. Adams         ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr);
11692e68589bSMark F. Adams       }
11702e68589bSMark F. Adams     }
11712e68589bSMark F. Adams     ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr);
11722e68589bSMark F. Adams   }
11732e68589bSMark F. Adams   else {
1174c8b0795cSMark F. Adams     nbnodes = bs*nloc;
1175c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11762e68589bSMark F. Adams   }
11772e68589bSMark F. Adams 
11782e68589bSMark F. Adams   /* get P0 */
11792e68589bSMark F. Adams   if( npe > 1 ){
11802e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
11812e68589bSMark F. Adams     PetscInt nnodes;
11822e68589bSMark F. Adams 
11832e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr);
11842e68589bSMark F. Adams     for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
11850cbbd2e1SMark F. Adams     ierr = PCGAMGGetDataWithGhosts( Gmat, 1, fid_glid_loc, &nnodes, &fiddata );
11862e68589bSMark F. Adams     CHKERRQ(ierr);
11872e68589bSMark F. Adams     ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
11882e68589bSMark F. Adams     for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11892e68589bSMark F. Adams     ierr = PetscFree( fiddata ); CHKERRQ(ierr);
1190c8b0795cSMark F. Adams     assert(nnodes==nbnodes/bs);
11912e68589bSMark F. Adams     ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr);
11922e68589bSMark F. Adams   }
11932e68589bSMark F. Adams   else {
11942e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
11952e68589bSMark F. Adams     for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk;
11962e68589bSMark F. Adams   }
11970cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11980cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11990cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
12002e68589bSMark F. Adams #endif
1201c8b0795cSMark F. Adams   {
1202c8b0795cSMark F. Adams     PetscReal *data_out;
12030cbbd2e1SMark F. Adams     ierr = formProl0( agg_lists, bs, data_cols, myCrs0, nbnodes,
1204c8b0795cSMark F. Adams                       data_w_ghost, flid_fgid, &data_out, Prol );
12052e68589bSMark F. Adams     CHKERRQ(ierr);
1206c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
1207c8b0795cSMark F. Adams     pc_gamg->data = data_out;
1208c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1209c8b0795cSMark F. Adams     pc_gamg->data_sz = data_cols*data_cols*nLocalSelected;
1210c8b0795cSMark F. Adams   }
12110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12120cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1213c8b0795cSMark F. Adams #endif
12142e68589bSMark F. Adams   if (npe > 1) ierr = PetscFree( data_w_ghost );      CHKERRQ(ierr);
12152e68589bSMark F. Adams   ierr = PetscFree( flid_fgid ); CHKERRQ(ierr);
12162e68589bSMark F. Adams 
1217c8b0795cSMark F. Adams   /* attach block size of columns */
1218c8b0795cSMark F. Adams   if( pc_gamg->col_bs_id == -1 ) {
1219c8b0795cSMark F. Adams     ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id ); assert(pc_gamg->col_bs_id != -1 );
1220c8b0795cSMark F. Adams   }
1221c8b0795cSMark F. Adams   ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr);
1222c8b0795cSMark F. Adams 
1223c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
12240cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12250cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
12260cbbd2e1SMark F. Adams #endif
1227c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1228c8b0795cSMark F. Adams }
1229c8b0795cSMark F. Adams 
1230c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1231c8b0795cSMark F. Adams /*
12320cbbd2e1SMark F. Adams    PCGAMGOptprol_AGG
1233c8b0795cSMark F. Adams 
1234c8b0795cSMark F. Adams   Input Parameter:
1235c8b0795cSMark F. Adams    . pc - this
1236c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1237c8b0795cSMark F. Adams  In/Output Parameter:
1238c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1239c8b0795cSMark F. Adams */
1240c8b0795cSMark F. Adams #undef __FUNCT__
12410cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG"
12420cbbd2e1SMark F. Adams PetscErrorCode PCGAMGOptprol_AGG( PC pc,
1243c8b0795cSMark F. Adams                                   const Mat Amat,
1244c8b0795cSMark F. Adams                                   Mat *a_P
1245c8b0795cSMark F. Adams                                   )
1246c8b0795cSMark F. Adams {
1247c8b0795cSMark F. Adams   PetscErrorCode ierr;
1248c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
1249c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1250c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
1251c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1252c8b0795cSMark F. Adams   PetscInt       jj;
1253c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
1254c8b0795cSMark F. Adams   Mat            Prol = *a_P;
1255c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1256c8b0795cSMark F. Adams 
1257c8b0795cSMark F. Adams   PetscFunctionBegin;
12580cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12590cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
12600cbbd2e1SMark F. Adams #endif
1261c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
1262c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
1263c8b0795cSMark F. Adams 
12642e68589bSMark F. Adams   /* smooth P0 */
1265c8b0795cSMark F. Adams   for( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){
12662e68589bSMark F. Adams     Mat tMat;
12672e68589bSMark F. Adams     Vec diag;
12682e68589bSMark F. Adams     PetscReal alpha, emax, emin;
12690cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12700cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12712e68589bSMark F. Adams #endif
12722e68589bSMark F. Adams     if( jj == 0 ) {
12732e68589bSMark F. Adams       KSP eksp;
12742e68589bSMark F. Adams       Vec bb, xx;
12752e68589bSMark F. Adams       PC pc;
12762e68589bSMark F. Adams       ierr = MatGetVecs( Amat, &bb, 0 );         CHKERRQ(ierr);
12772e68589bSMark F. Adams       ierr = MatGetVecs( Amat, &xx, 0 );         CHKERRQ(ierr);
12782e68589bSMark F. Adams       {
12792e68589bSMark F. Adams         PetscRandom    rctx;
12802e68589bSMark F. Adams         ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
12812e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
12822e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
12832e68589bSMark F. Adams         ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
12842e68589bSMark F. Adams       }
12852e68589bSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);                            CHKERRQ(ierr);
12862e68589bSMark F. Adams       ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
12872e68589bSMark F. Adams       ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
12882e68589bSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE );    CHKERRQ(ierr);
12892e68589bSMark F. Adams       ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN );
12902e68589bSMark F. Adams       CHKERRQ( ierr );
12912e68589bSMark F. Adams       ierr = KSPGetPC( eksp, &pc );                              CHKERRQ( ierr );
12922e68589bSMark F. Adams       ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr);  /* smoother */
12932e68589bSMark F. Adams       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);
12942e68589bSMark F. Adams       CHKERRQ(ierr);
12952e68589bSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
12962e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE );        CHKERRQ(ierr);
12972e68589bSMark F. Adams 
12982e68589bSMark F. Adams       /* solve - keep stuff out of logging */
12992e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
13002e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
13012e68589bSMark F. Adams       ierr = KSPSolve( eksp, bb, xx );                              CHKERRQ(ierr);
13022e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
13032e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
13042e68589bSMark F. Adams 
13052e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
13062e68589bSMark F. Adams       if( verbose ) {
1307c8b0795cSMark F. Adams         PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
13082e68589bSMark F. Adams                     __FUNCT__,emax,emin,PCJACOBI);
13092e68589bSMark F. Adams       }
13102e68589bSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
13112e68589bSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
13122e68589bSMark F. Adams       ierr = KSPDestroy( &eksp );     CHKERRQ(ierr);
13132e68589bSMark F. Adams 
13142e68589bSMark F. Adams       if( pc_gamg->emax_id == -1 ) {
13152e68589bSMark F. Adams         ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id );
13162e68589bSMark F. Adams         assert(pc_gamg->emax_id != -1 );
13172e68589bSMark F. Adams       }
13182e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr);
13192e68589bSMark F. Adams     }
13202e68589bSMark F. Adams 
13212e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
13222e68589bSMark F. Adams     ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat );   CHKERRQ(ierr);
13232e68589bSMark F. Adams     ierr = MatGetVecs( Amat, &diag, 0 );    CHKERRQ(ierr);
13242e68589bSMark F. Adams     ierr = MatGetDiagonal( Amat, diag );    CHKERRQ(ierr); /* effectively PCJACOBI */
13252e68589bSMark F. Adams     ierr = VecReciprocal( diag );         CHKERRQ(ierr);
13262e68589bSMark F. Adams     ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr);
13272e68589bSMark F. Adams     ierr = VecDestroy( &diag );           CHKERRQ(ierr);
13282e68589bSMark F. Adams     alpha = -1.5/emax;
13292e68589bSMark F. Adams     ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN );           CHKERRQ(ierr);
13302e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
13312e68589bSMark F. Adams     Prol = tMat;
13320cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13330cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
13342e68589bSMark F. Adams #endif
13352e68589bSMark F. Adams   }
13360cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13370cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
13380cbbd2e1SMark F. Adams #endif
1339c8b0795cSMark F. Adams   *a_P = Prol;
1340c8b0795cSMark F. Adams 
1341c8b0795cSMark F. Adams   PetscFunctionReturn(0);
13422e68589bSMark F. Adams }
13432e68589bSMark F. Adams 
1344c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1345c8b0795cSMark F. Adams /*
1346c8b0795cSMark F. Adams    PCCreateGAMG_AGG
13472e68589bSMark F. Adams 
1348c8b0795cSMark F. Adams   Input Parameter:
1349c8b0795cSMark F. Adams    . pc -
1350c8b0795cSMark F. Adams */
1351c8b0795cSMark F. Adams #undef __FUNCT__
1352c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1353c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG( PC pc )
1354c8b0795cSMark F. Adams {
1355c8b0795cSMark F. Adams   PetscErrorCode  ierr;
1356c8b0795cSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1357c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1358c8b0795cSMark F. Adams   PC_GAMG_AGG      *pc_gamg_agg;
13592e68589bSMark F. Adams 
1360c8b0795cSMark F. Adams   PetscFunctionBegin;
1361c8b0795cSMark F. Adams   /* create sub context for SA */
1362c8b0795cSMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr);
1363c8b0795cSMark F. Adams   assert(!pc_gamg->subctx);
1364c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1365c8b0795cSMark F. Adams 
1366c8b0795cSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1367c8b0795cSMark F. Adams   pc->ops->destroy        = PCDestroy_AGG;
1368c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1369c8b0795cSMark F. Adams 
1370c8b0795cSMark F. Adams   /* set internal function pointers */
1371c8b0795cSMark F. Adams   pc_gamg->graph = PCGAMGgraph_AGG;
1372b43b03e9SMark F. Adams   pc_gamg->coarsen = PCGAMGCoarsen_AGG;
13730cbbd2e1SMark F. Adams   pc_gamg->prolongator = PCGAMGProlongator_AGG;
13740cbbd2e1SMark F. Adams   pc_gamg->optprol = PCGAMGOptprol_AGG;
1375c8b0795cSMark F. Adams 
1376c8b0795cSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_AGG;
1377c8b0795cSMark F. Adams 
1378c8b0795cSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1379c8b0795cSMark F. Adams                                             "PCSetCoordinates_C",
1380c8b0795cSMark F. Adams                                             "PCSetCoordinates_AGG",
1381c8b0795cSMark F. Adams                                             PCSetCoordinates_AGG);
13822e68589bSMark F. Adams   PetscFunctionReturn(0);
13832e68589bSMark F. Adams }
1384