xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision e696ed0bff0adf00a24c551cb9f31e65f9331d7c)
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
238302f38e8SMark F. Adams      - collective
2392e68589bSMark F. Adams 
2402e68589bSMark F. Adams    Input Parameter:
2412e68589bSMark F. Adams    . pc - the preconditioner context
242a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
243302f38e8SMark F. Adams    . a_nloc - number of vertices local
244302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2452e68589bSMark F. Adams */
2462e68589bSMark F. Adams EXTERN_C_BEGIN
2472e68589bSMark F. Adams #undef __FUNCT__
2482e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
249302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords )
2502e68589bSMark F. Adams {
2512e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
2522e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2532e68589bSMark F. Adams   PetscErrorCode  ierr;
25469344418SMark F. Adams   PetscInt        arrsz,kk,ii,jj,nloc,ndatarows,ndf;
255a2f3521dSMark F. Adams   Mat             mat = pc->pmat;
256a2f3521dSMark F. Adams   /* MPI_Comm     wcomm = ((PetscObject)pc)->comm; */
2572e68589bSMark F. Adams 
2582e68589bSMark F. Adams   PetscFunctionBegin;
259a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
260a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
261302f38e8SMark F. Adams   nloc = a_nloc;
2622e68589bSMark F. Adams 
2632e68589bSMark F. Adams   /* SA: null space vectors */
26469344418SMark F. Adams   ierr = MatGetBlockSize( mat, &ndf ); CHKERRQ( ierr ); /* this does not work for Stokes */
26569344418SMark F. Adams   if ( coords && ndf==1 ) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
266a2f3521dSMark F. Adams   else if ( coords ) {
267c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
26869344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
26969344418SMark F. Adams     if (ndm != ndf) {
270c666adbfSMark 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);
271a2f3521dSMark F. Adams     }
27269344418SMark F. Adams   }
27369344418SMark F. Adams   else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
27469344418SMark F. Adams   pc_gamg->data_cell_rows = ndatarows = ndf;     assert(pc_gamg->data_cell_cols>0);
275c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2762e68589bSMark F. Adams 
2772e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2782e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
2792e68589bSMark F. Adams     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
280302f38e8SMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
2812e68589bSMark F. Adams   }
2822e68589bSMark F. Adams   /* copy data in - column oriented */
2832e68589bSMark F. Adams   for (kk=0;kk<nloc;kk++){
28469344418SMark F. Adams     const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */
28569344418SMark F. Adams     PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
286c8b0795cSMark F. Adams     if ( pc_gamg->data_cell_cols==1 ) *data = 1.0;
2872e68589bSMark F. Adams     else {
28869344418SMark F. Adams       /* translational modes */
289a2f3521dSMark F. Adams       for (ii=0;ii<ndatarows;ii++)
290a2f3521dSMark F. Adams         for (jj=0;jj<ndatarows;jj++)
29169344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2922e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
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];
3002e68589bSMark F. Adams         }
3012e68589bSMark F. Adams         else {
3022e68589bSMark F. Adams           data += 3*M;
3032e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
3042e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
3052e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
3062e68589bSMark F. Adams         }
3072e68589bSMark F. Adams       }
3082e68589bSMark F. Adams     }
3092e68589bSMark F. Adams   }
3102e68589bSMark F. Adams 
3112e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
3122e68589bSMark F. Adams 
3132e68589bSMark F. Adams   PetscFunctionReturn(0);
3142e68589bSMark F. Adams }
3152e68589bSMark F. Adams EXTERN_C_END
3162e68589bSMark F. Adams 
317b43b03e9SMark F. Adams typedef PetscInt NState;
318b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
319b43b03e9SMark F. Adams static const NState DELETED=-1;
320b43b03e9SMark F. Adams static const NState REMOVED=-3;
321b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
322b43b03e9SMark F. Adams 
323c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
324c8b0795cSMark F. Adams /*
325b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
326b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
327c8b0795cSMark F. Adams 
328c8b0795cSMark F. Adams    Input Parameter:
329c8b0795cSMark F. Adams    . Gmat_2 - glabal matrix of graph (data not defined)
330c8b0795cSMark F. Adams    . Gmat_1 - base graph to grab with
331c8b0795cSMark F. Adams    Input/Output Parameter:
3320cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids )
333c8b0795cSMark F. Adams */
334c8b0795cSMark F. Adams #undef __FUNCT__
335c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3360cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */
3370cbbd2e1SMark F. Adams                                   const Mat Gmat_1, /* base graph */
3380cbbd2e1SMark F. Adams                                   /* const IS selected_2, [nselected local] selected vertices */
3390cbbd2e1SMark F. Adams                                   PetscCoarsenData *aggs_2 /* [nselected local] global ID of aggregate */
340c8b0795cSMark F. Adams                                   )
341c8b0795cSMark F. Adams {
342c8b0795cSMark F. Adams   PetscErrorCode ierr;
343c8b0795cSMark F. Adams   PetscBool      isMPI;
344c8b0795cSMark F. Adams   Mat_SeqAIJ    *matA_1, *matB_1=0, *matA_2, *matB_2=0;
345c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat_2)->comm;
3460cbbd2e1SMark F. Adams   PetscMPIInt    mype,npe;
3470cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
348c8b0795cSMark F. Adams   Mat_MPIAIJ    *mpimat_2 = 0, *mpimat_1=0;
349c8b0795cSMark F. Adams   const PetscInt nloc = Gmat_2->rmap->n;
3500cbbd2e1SMark F. Adams   PetscScalar   *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3510cbbd2e1SMark F. Adams   PetscInt      *lid_cprowID_1;
352c8b0795cSMark F. Adams   NState        *lid_state;
3530cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
354c8b0795cSMark F. Adams 
355c8b0795cSMark F. Adams   PetscFunctionBegin;
356c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype );   CHKERRQ(ierr);
3570cbbd2e1SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe );   CHKERRQ(ierr);
358c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);  CHKERRQ(ierr);
359c8b0795cSMark F. Adams 
3600cbbd2e1SMark F. Adams   if ( PETSC_FALSE ) {
361c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
362c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
363c8b0795cSMark F. Adams     PetscViewerASCIIOpen(wcomm,fname,&viewer);
364c8b0795cSMark F. Adams     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
365c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr);
366c8b0795cSMark F. Adams     ierr = PetscViewerDestroy( &viewer );
367c8b0795cSMark F. Adams   }
368c8b0795cSMark F. Adams 
369c8b0795cSMark F. Adams   /* get submatrices */
370251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr);
371c8b0795cSMark F. Adams   if (isMPI) {
372c8b0795cSMark F. Adams     /* grab matrix objects */
373c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
374c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
375c8b0795cSMark F. Adams     matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data;
376c8b0795cSMark F. Adams     matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data;
377c8b0795cSMark F. Adams     matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data;
378c8b0795cSMark F. Adams     matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data;
379c8b0795cSMark F. Adams 
380c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
381c8b0795cSMark F. Adams     matB_1->compressedrow.check = PETSC_TRUE;
382c8b0795cSMark F. Adams     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);
383c8b0795cSMark F. Adams     CHKERRQ(ierr);
384c8b0795cSMark F. Adams 
385c8b0795cSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr);
3860cbbd2e1SMark F. Adams     for ( lid = 0 ; lid < nloc ; lid++ ) lid_cprowID_1[lid] = -1;
387c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
388c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
389c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
390c8b0795cSMark F. Adams     }
391c8b0795cSMark F. Adams   }
392c8b0795cSMark F. Adams   else {
393c8b0795cSMark F. Adams     matA_1 = (Mat_SeqAIJ*)Gmat_1->data;
394c8b0795cSMark F. Adams     matA_2 = (Mat_SeqAIJ*)Gmat_2->data;
3950cbbd2e1SMark F. Adams     lid_cprowID_1 = PETSC_NULL;
396c8b0795cSMark F. Adams   }
397c8b0795cSMark F. Adams   assert( matA_1 && !matA_1->compressedrow.use );
398c8b0795cSMark F. Adams   assert( matB_1==0 || matB_1->compressedrow.use );
399c8b0795cSMark F. Adams   assert( matA_2 && !matA_2->compressedrow.use );
400c8b0795cSMark F. Adams   assert( matB_2==0 || matB_2->compressedrow.use );
401c8b0795cSMark F. Adams 
402c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
403c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr);
4040cbbd2e1SMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscScalar), &lid_parent_gid ); CHKERRQ(ierr);
405c8b0795cSMark F. Adams   for ( lid = 0 ; lid < nloc ; lid++ ) {
4060cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
407c8b0795cSMark F. Adams     lid_state[lid] = DELETED;
408c8b0795cSMark F. Adams   }
4090cbbd2e1SMark F. Adams 
4100cbbd2e1SMark F. Adams   /* set lid_state */
4110cbbd2e1SMark F. Adams   for ( lid = 0 ; lid < nloc ; lid++ ) {
41241b27cdeSMark F. Adams     PetscCDPos pos;
413e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
414e78576d6SMark F. Adams     if ( pos ) {
415e78576d6SMark F. Adams       PetscInt gid1;
416ffc955d6SMark F. Adams       ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr); assert(gid1==lid+my0);
4170cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
418b43b03e9SMark F. Adams     }
419b43b03e9SMark F. Adams   }
4200cbbd2e1SMark F. Adams 
4210cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
422c8b0795cSMark F. Adams   for (lid=kk=0;lid<nloc;lid++){
423c8b0795cSMark F. Adams     NState state = lid_state[lid];
424c8b0795cSMark F. Adams     if ( IS_SELECTED(state) ){
42541b27cdeSMark F. Adams       PetscCDPos pos;
426e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
427e78576d6SMark F. Adams       while(pos){
428e78576d6SMark F. Adams         PetscInt gid1;
429ffc955d6SMark F. Adams         ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr);
430e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr);
431e78576d6SMark F. Adams 
4320cbbd2e1SMark F. Adams         if ( gid1 >= my0 && gid1 < Iend ){
4330cbbd2e1SMark F. Adams           lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
434c8b0795cSMark F. Adams         }
435c8b0795cSMark F. Adams       }
4360cbbd2e1SMark F. Adams     }
4370cbbd2e1SMark F. Adams   }
4380cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
439c8b0795cSMark F. Adams   if (isMPI) {
440c8b0795cSMark F. Adams     Vec          tempVec;
441c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
442c8b0795cSMark F. Adams     ierr = MatGetVecs( Gmat_1, &tempVec, 0 );         CHKERRQ(ierr);
443c8b0795cSMark F. Adams     for (kk=0,j=my0;kk<nloc;kk++,j++){
444c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
445c8b0795cSMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
446c8b0795cSMark F. Adams     }
447c8b0795cSMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
448c8b0795cSMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
449c8b0795cSMark F. Adams     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
450c8b0795cSMark F. Adams     CHKERRQ(ierr);
451c8b0795cSMark F. Adams     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
452c8b0795cSMark F. Adams     CHKERRQ(ierr);
453c8b0795cSMark F. Adams     ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
454c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
455c8b0795cSMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
456c8b0795cSMark F. Adams     CHKERRQ(ierr);
457c8b0795cSMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
458c8b0795cSMark F. Adams     CHKERRQ(ierr);
459c8b0795cSMark F. Adams     ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
4600cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4610cbbd2e1SMark F. Adams     for (kk=0,j=my0;kk<nloc;kk++,j++){
4620cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4630cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
4640cbbd2e1SMark F. Adams     }
4650cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
4660cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
4670cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghost_par_orig2 ); CHKERRQ(ierr);
4680cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);
4690cbbd2e1SMark F. Adams     CHKERRQ(ierr);
4700cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);
4710cbbd2e1SMark F. Adams     CHKERRQ(ierr);
4720cbbd2e1SMark F. Adams     ierr = VecGetArray( ghost_par_orig2, &cpcol_2_par_orig ); CHKERRQ(ierr);
4730cbbd2e1SMark F. Adams 
474c8b0795cSMark F. Adams     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
475c8b0795cSMark F. Adams   } /* ismpi */
476c8b0795cSMark F. Adams 
477c8b0795cSMark F. Adams   /* doit */
478c8b0795cSMark F. Adams   for (lid=0;lid<nloc;lid++){
479c8b0795cSMark F. Adams     NState state = lid_state[lid];
4800cbbd2e1SMark F. Adams     if ( IS_SELECTED(state) ) {
4810cbbd2e1SMark F. Adams       /* steal locals */
482c8b0795cSMark F. Adams       ii = matA_1->i; n = ii[lid+1] - ii[lid];
483c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
484c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4850cbbd2e1SMark F. Adams         PetscInt lidj = idx[j], sgid;
486c8b0795cSMark F. Adams         NState statej = lid_state[lidj];
4870cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4880cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4890cbbd2e1SMark F. Adams           if ( sgid >= my0 && sgid < Iend ){       /* I'm stealing this local from a local sgid */
4900cbbd2e1SMark F. Adams             PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0;
49141b27cdeSMark F. Adams             PetscCDPos pos,last=PETSC_NULL;
492c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
493e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos); CHKERRQ(ierr);
494e78576d6SMark F. Adams             while(pos){
495e78576d6SMark F. Adams               PetscInt gid;
496ffc955d6SMark F. Adams               ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr);
4970cbbd2e1SMark F. Adams               if ( gid == gidj ) {
4980cbbd2e1SMark F. Adams                 assert(last);
49941b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode( aggs_2, slid, last ); CHKERRQ(ierr);
50041b27cdeSMark F. Adams                 ierr = PetscCDAppendNode( aggs_2, lid, pos );       CHKERRQ(ierr);
5010cbbd2e1SMark F. Adams                 hav = 1;
502c8b0795cSMark F. Adams                 break;
503c8b0795cSMark F. Adams               }
5040cbbd2e1SMark F. Adams               else last = pos;
505e78576d6SMark F. Adams 
506e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos); CHKERRQ(ierr);
507c8b0795cSMark F. Adams             }
508c8b0795cSMark F. Adams             if (hav!=1){
509c666adbfSMark F. Adams               if (hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
510c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
511c8b0795cSMark F. Adams             }
512c8b0795cSMark F. Adams           }
5130cbbd2e1SMark F. Adams           else{            /* I'm stealing this local, owned by a ghost */
514c666adbfSMark F. Adams             if (sgid != -1 ) {
515c666adbfSMark F. Adams               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.");
516c666adbfSMark F. Adams             }
51741b27cdeSMark F. Adams             ierr = PetscCDAppendID( aggs_2, lid, lidj+my0 );      CHKERRQ(ierr);
518c8b0795cSMark F. Adams           }
519c8b0795cSMark F. Adams         }
5200cbbd2e1SMark F. Adams       } /* local neighbors */
521c8b0795cSMark F. Adams     }
522c8b0795cSMark F. Adams     else if ( state == DELETED && lid_cprowID_1 ) {
5230cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
524c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
525c8b0795cSMark F. Adams       if ( (ix=lid_cprowID_1[lid]) != -1 ){
526c8b0795cSMark F. Adams         ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
527c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
528c8b0795cSMark F. Adams         for ( j=0 ; j<n ; j++ ) {
529c8b0795cSMark F. Adams           PetscInt cpid = idx[j];
530c8b0795cSMark F. Adams           NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
531c8b0795cSMark F. Adams           if ( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */
5320cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
5330cbbd2e1SMark F. Adams             if ( sgidold>=my0 && sgidold<Iend ) { /* this was mine */
5340cbbd2e1SMark F. Adams               PetscInt hav=0,oldslidj=sgidold-my0;
53541b27cdeSMark F. Adams               PetscCDPos pos,last=PETSC_NULL;
5360cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
537e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr);
538e78576d6SMark F. Adams               while( pos ) {
539e78576d6SMark F. Adams                 PetscInt gid;
540ffc955d6SMark F. Adams                 ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr);
5410cbbd2e1SMark F. Adams                 if ( lid+my0 == gid ) {
5420cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
5430cbbd2e1SMark F. Adams                   assert(last);
54441b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode( aggs_2, oldslidj, last ); CHKERRQ(ierr);
5450cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
5460cbbd2e1SMark F. Adams                   hav = 1;
547c8b0795cSMark F. Adams                   break;
548c8b0795cSMark F. Adams                 }
5490cbbd2e1SMark F. Adams                 else last = pos;
550e78576d6SMark F. Adams 
551e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr);
552c8b0795cSMark F. Adams               }
553c8b0795cSMark F. Adams               if (hav!=1){
554c666adbfSMark F. Adams                 if (hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
555c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
556c8b0795cSMark F. Adams               }
557c8b0795cSMark F. Adams             }
5580cbbd2e1SMark F. Adams             else {
5590cbbd2e1SMark F. Adams               /* ghosts remove this later */
5600cbbd2e1SMark F. Adams             }
561c8b0795cSMark F. Adams           }
562c8b0795cSMark F. Adams         }
563c8b0795cSMark F. Adams       }
564c8b0795cSMark F. Adams     } /* selected/deleted */
565c8b0795cSMark F. Adams   } /* node loop */
566c8b0795cSMark F. Adams 
567c8b0795cSMark F. Adams   if ( isMPI ) {
5680cbbd2e1SMark F. Adams     PetscScalar *cpcol_2_parent,*cpcol_2_gid;
5690cbbd2e1SMark F. Adams     Vec          tempVec,ghostgids2,ghostparents2;
5700cbbd2e1SMark F. Adams     PetscInt     cpid,nghost_2;
5710cbbd2e1SMark F. Adams     GAMGHashTable gid_cpid;
572c8b0795cSMark F. Adams 
5730cbbd2e1SMark F. Adams     ierr = VecGetSize( mpimat_2->lvec, &nghost_2 );   CHKERRQ(ierr);
574c8b0795cSMark F. Adams     ierr = MatGetVecs( Gmat_2, &tempVec, 0 );         CHKERRQ(ierr);
5750cbbd2e1SMark F. Adams 
5760cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
577c8b0795cSMark F. Adams     for (kk=0,j=my0;kk<nloc;kk++,j++){
5780cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES );  CHKERRQ(ierr);
579c8b0795cSMark F. Adams     }
580c8b0795cSMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
581c8b0795cSMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
5820cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghostparents2 ); CHKERRQ(ierr);
5830cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);
584c8b0795cSMark F. Adams     CHKERRQ(ierr);
5850cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);
586c8b0795cSMark F. Adams     CHKERRQ(ierr);
5870cbbd2e1SMark F. Adams     ierr = VecGetArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr);
5880cbbd2e1SMark F. Adams 
5890cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5900cbbd2e1SMark F. Adams     for (kk=0,j=my0;kk<nloc;kk++,j++){
5910cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5920cbbd2e1SMark F. Adams       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
5930cbbd2e1SMark F. Adams     }
5940cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
5950cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
5960cbbd2e1SMark F. Adams     ierr = VecDuplicate( mpimat_2->lvec, &ghostgids2 ); CHKERRQ(ierr);
5970cbbd2e1SMark F. Adams     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);
5980cbbd2e1SMark F. Adams     CHKERRQ(ierr);
5990cbbd2e1SMark F. Adams     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);
6000cbbd2e1SMark F. Adams     CHKERRQ(ierr);
6010cbbd2e1SMark F. Adams     ierr = VecGetArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr);
6020cbbd2e1SMark F. Adams 
603c8b0795cSMark F. Adams     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
604c8b0795cSMark F. Adams 
6050cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
6060cbbd2e1SMark F. Adams     ierr = GAMGTableCreate( 2*nghost_2, &gid_cpid ); CHKERRQ(ierr);
6070cbbd2e1SMark F. Adams     for ( cpid = 0 ; cpid < nghost_2 ; cpid++ ) {
6080cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
6090cbbd2e1SMark F. Adams       if ( state==DELETED ) {
6100cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6110cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
6120cbbd2e1SMark F. Adams         if ( sgid_old == -1 && sgid_new != -1 ) {
6130cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6140cbbd2e1SMark F. Adams           ierr = GAMGTableAdd( &gid_cpid, gid, cpid ); CHKERRQ(ierr);
6150cbbd2e1SMark F. Adams         }
6160cbbd2e1SMark F. Adams       }
6170cbbd2e1SMark F. Adams     }
618c8b0795cSMark F. Adams 
6190cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
620c8b0795cSMark F. Adams     for (lid=0;lid<nloc;lid++){
621c8b0795cSMark F. Adams       NState state = lid_state[lid];
622c8b0795cSMark F. Adams       if ( IS_SELECTED(state) ){
62341b27cdeSMark F. Adams         PetscCDPos pos,last=PETSC_NULL;
624c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
625e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr);
626e78576d6SMark F. Adams         while(pos){
627e78576d6SMark F. Adams           PetscInt gid;
628ffc955d6SMark F. Adams           ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr);
629e78576d6SMark F. Adams 
6300cbbd2e1SMark F. Adams           if ( gid < my0 || gid >= Iend ) {
6310cbbd2e1SMark F. Adams             ierr = GAMGTableFind( &gid_cpid, gid, &cpid ); CHKERRQ(ierr);
6320cbbd2e1SMark F. Adams             if ( cpid != -1 ) {
6330cbbd2e1SMark F. Adams               /* a moved ghost - */
6340cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
63541b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode( aggs_2, lid, last ); CHKERRQ(ierr);
6360cbbd2e1SMark F. Adams             }
6370cbbd2e1SMark F. Adams             else last = pos;
6380cbbd2e1SMark F. Adams           }
6390cbbd2e1SMark F. Adams           else last = pos;
640e78576d6SMark F. Adams 
641e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr);
642c8b0795cSMark F. Adams         } /* loop over list of deleted */
643c8b0795cSMark F. Adams       } /* selected */
644c8b0795cSMark F. Adams     }
6450cbbd2e1SMark F. Adams     ierr = GAMGTableDestroy( &gid_cpid ); CHKERRQ(ierr);
646c8b0795cSMark F. Adams 
6470cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
6480cbbd2e1SMark F. Adams     for ( cpid = 0 ; cpid < nghost_2 ; cpid++ ){
6490cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
6500cbbd2e1SMark F. Adams       if ( sgid_new >= my0 && sgid_new < Iend ) { /* this is mine */
6510cbbd2e1SMark F. Adams         PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
6520cbbd2e1SMark F. Adams         PetscInt slid_new=sgid_new-my0,hav=0;
65341b27cdeSMark F. Adams         PetscCDPos pos;
6540cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
655e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos); CHKERRQ(ierr);
656e78576d6SMark F. Adams         while(pos){
657e78576d6SMark F. Adams           PetscInt gidj;
658ffc955d6SMark F. Adams           ierr = PetscLLNGetID( pos, &gidj ); CHKERRQ(ierr);
659e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos); CHKERRQ(ierr);
660e78576d6SMark F. Adams 
6610cbbd2e1SMark F. Adams           if ( gidj == gid ) { hav = 1; break; }
662c8b0795cSMark F. Adams         }
663c8b0795cSMark F. Adams         if ( hav != 1 ){
664ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
66541b27cdeSMark F. Adams           ierr = PetscCDAppendID( aggs_2, slid_new, gid );      CHKERRQ(ierr);
666c8b0795cSMark F. Adams         }
667c8b0795cSMark F. Adams       }
668c8b0795cSMark F. Adams     }
669c8b0795cSMark F. Adams 
6700cbbd2e1SMark F. Adams     ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
6710cbbd2e1SMark F. Adams     ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
6720cbbd2e1SMark F. Adams     ierr = VecRestoreArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr);
6730cbbd2e1SMark F. Adams     ierr = VecRestoreArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr);
674c8b0795cSMark F. Adams     ierr = PetscFree( lid_cprowID_1 );  CHKERRQ(ierr);
6750cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghostgids2 ); CHKERRQ(ierr);
6760cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghostparents2 ); CHKERRQ(ierr);
6770cbbd2e1SMark F. Adams     ierr = VecDestroy( &ghost_par_orig2 ); CHKERRQ(ierr);
678c8b0795cSMark F. Adams   }
679c8b0795cSMark F. Adams 
6800cbbd2e1SMark F. Adams   ierr = PetscFree( lid_parent_gid );  CHKERRQ(ierr);
681c8b0795cSMark F. Adams   ierr = PetscFree( lid_state );  CHKERRQ(ierr);
682c8b0795cSMark F. Adams 
683c8b0795cSMark F. Adams   PetscFunctionReturn(0);
684c8b0795cSMark F. Adams }
6852e68589bSMark F. Adams 
6862e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6872e68589bSMark F. Adams /*
688a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
689a2f3521dSMark F. Adams       Looks in Mat for near null space.
690a2f3521dSMark F. Adams       Does not work for Stokes
6912e68589bSMark F. Adams 
6922e68589bSMark F. Adams   Input Parameter:
6932e68589bSMark F. Adams    . pc -
694a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6952e68589bSMark F. Adams */
6962e68589bSMark F. Adams #undef __FUNCT__
6972e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
698b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG( PC pc, Mat a_A )
6992e68589bSMark F. Adams {
7002e68589bSMark F. Adams   PetscErrorCode  ierr;
701b8cd405aSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
702b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
703b8cd405aSMark F. Adams   MatNullSpace mnull;
704b8cd405aSMark F. Adams 
7052e68589bSMark F. Adams   PetscFunctionBegin;
706b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace( a_A, &mnull ); CHKERRQ(ierr);
707b8cd405aSMark F. Adams   if ( !mnull ) {
708a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
7099d1b15c3SMark F. Adams     ierr = MatGetBlockSize( a_A, &bs ); CHKERRQ( ierr );
7109d1b15c3SMark F. Adams     ierr = MatGetLocalSize( a_A, &MM, &NN ); CHKERRQ( ierr );  assert(MM%bs==0);
711a2f3521dSMark F. Adams     ierr = PCSetCoordinates_AGG( pc, bs, MM/bs, PETSC_NULL ); CHKERRQ(ierr);
712b8cd405aSMark F. Adams   }
713b8cd405aSMark F. Adams   else {
714b8cd405aSMark F. Adams     PetscReal *nullvec;
715b8cd405aSMark F. Adams     PetscBool has_const;
716b8cd405aSMark F. Adams     PetscInt i,j,mlocal,nvec,bs;
717b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
718b8cd405aSMark F. Adams     ierr = MatGetLocalSize(a_A,&mlocal,PETSC_NULL);CHKERRQ(ierr);
719b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs( mnull, &has_const, &nvec, &vecs ); CHKERRQ(ierr);
7208caf3d72SBarry Smith     ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
721b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
722b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
723b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
724b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
725b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
726b8cd405aSMark F. Adams     }
727b8cd405aSMark F. Adams     pc_gamg->data = nullvec;
728b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
729a2f3521dSMark F. Adams     ierr = MatGetBlockSize( a_A, &bs ); CHKERRQ( ierr ); /* this does not work for Stokes */
730b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
731b8cd405aSMark F. Adams   }
7322e68589bSMark F. Adams   PetscFunctionReturn(0);
7332e68589bSMark F. Adams }
7342e68589bSMark F. Adams 
7352e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
7362e68589bSMark F. Adams /*
7372e68589bSMark F. Adams  formProl0
7382e68589bSMark F. Adams 
7392e68589bSMark F. Adams    Input Parameter:
7400cbbd2e1SMark F. Adams    . agg_llists - list of arrays with aggregates
7412e68589bSMark F. Adams    . bs - block size
7420cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7430cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7442e68589bSMark F. Adams    . data_stride - bs*(nloc nodes + ghost nodes)
7452e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7462e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7472e68589bSMark F. Adams   Output Parameter:
7482e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7492e68589bSMark F. Adams    . a_Prol - prolongation operator
7502e68589bSMark F. Adams */
7512e68589bSMark F. Adams #undef __FUNCT__
7522e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7530cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists,/* list from selected vertices of aggregate unselected vertices */
7540cbbd2e1SMark F. Adams                                 const PetscInt bs,          /* (row) block size */
7550cbbd2e1SMark F. Adams                                 const PetscInt nSAvec,      /* column bs */
7560cbbd2e1SMark F. Adams                                 const PetscInt my0crs,      /* global index of start of locals */
7570cbbd2e1SMark F. Adams                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
7580cbbd2e1SMark F. Adams                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
7590cbbd2e1SMark F. Adams                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
7602e68589bSMark F. Adams                                 PetscReal **a_data_out,
7612e68589bSMark F. Adams                                 Mat a_Prol /* prolongation operator (output)*/
7622e68589bSMark F. Adams                                 )
7632e68589bSMark F. Adams {
7642e68589bSMark F. Adams   PetscErrorCode ierr;
7650cbbd2e1SMark F. Adams   PetscInt  Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7662e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
7672e68589bSMark F. Adams   PetscMPIInt    mype, npe;
7682e68589bSMark F. Adams   PetscReal      *out_data;
76941b27cdeSMark F. Adams   PetscCDPos         pos;
7700cbbd2e1SMark F. Adams   GAMGHashTable  fgid_flid;
7710cbbd2e1SMark F. Adams 
772797e13b7SMark F. Adams /* #define OUT_AGGS */
7739057884aSMark F. Adams #ifdef OUT_AGGS
774f7620de1SMatthew G Knepley   static PetscInt llev = 0; char fname[32]; FILE *file = PETSC_NULL; PetscInt pM;
7759057884aSMark F. Adams #endif
7762e68589bSMark F. Adams 
7772e68589bSMark F. Adams   PetscFunctionBegin;
7782e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
7792e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
7802e68589bSMark F. Adams   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );    CHKERRQ(ierr);
7810cbbd2e1SMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
7820cbbd2e1SMark F. Adams   Iend /= bs;
7830cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7842e68589bSMark F. Adams 
7850cbbd2e1SMark F. Adams   ierr = GAMGTableCreate( 2*nghosts, &fgid_flid ); CHKERRQ(ierr);
7860cbbd2e1SMark F. Adams   for (kk=0;kk<nghosts;kk++) {
7870cbbd2e1SMark F. Adams     ierr = GAMGTableAdd( &fgid_flid, flid_fgid[nloc+kk], nloc+kk ); CHKERRQ(ierr);
7882e68589bSMark F. Adams   }
7892e68589bSMark F. Adams 
7900cbbd2e1SMark F. Adams #ifdef OUT_AGGS
7910cbbd2e1SMark F. Adams   sprintf(fname,"aggs_%d_%d.m",llev++,mype);
7920cbbd2e1SMark F. Adams   if (llev==1) {
7930cbbd2e1SMark F. Adams     file = fopen(fname,"w");
7940cbbd2e1SMark F. Adams   }
7950cbbd2e1SMark F. Adams   MatGetSize( a_Prol, &pM, &jj );
7960cbbd2e1SMark F. Adams #endif
7970cbbd2e1SMark F. Adams 
7980cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7990cbbd2e1SMark F. Adams   for (nSelected=mm=0;mm<nloc;mm++) {
800e78576d6SMark F. Adams     PetscBool ise;
801e78576d6SMark F. Adams     ierr = PetscCDEmptyAt( agg_llists, mm, &ise ); CHKERRQ(ierr);
802e78576d6SMark F. Adams     if ( !ise ) nSelected++;
8030cbbd2e1SMark F. Adams   }
8040cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn( a_Prol, &ii, &jj ); CHKERRQ(ierr);
8050cbbd2e1SMark F. Adams   assert((ii/nSAvec)==my0crs); assert(nSelected==(jj-ii)/nSAvec);
8060cbbd2e1SMark F. Adams 
8072e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
8080cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
8090cbbd2e1SMark F. Adams   ierr = PetscMalloc( out_data_stride*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr);
8100cbbd2e1SMark F. Adams   for (ii=0;ii<out_data_stride*nSAvec;ii++) {
8110cbbd2e1SMark F. Adams     out_data[ii]=1.e300;
8120cbbd2e1SMark F. Adams   }
8130cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
8142e68589bSMark F. Adams 
8152e68589bSMark F. Adams   /* find points and set prolongation */
816c8b0795cSMark F. Adams   minsz = 100;
8172e68589bSMark F. Adams   ndone = 0;
8180cbbd2e1SMark F. Adams   for ( mm = clid = 0 ; mm < nloc ; mm++ ){
819e78576d6SMark F. Adams     ierr = PetscCDSizeAt( agg_llists, mm, &jj ); CHKERRQ(ierr);
820e78576d6SMark F. Adams     if ( jj > 0 ) {
8210cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
8220cbbd2e1SMark F. Adams       PetscInt cids[100]; /* max bs */
8230cbbd2e1SMark F. Adams       PetscBLASInt asz=jj,M=asz*bs,N=nSAvec,INFO;
8242e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs;
8252e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
8262e68589bSMark F. Adams       PetscInt       *fids;
82765d7b583SSatish Balay       PetscReal      *data;
8280cbbd2e1SMark F. Adams       /* count agg */
8290cbbd2e1SMark F. Adams       if ( asz<minsz ) minsz = asz;
8300cbbd2e1SMark F. Adams 
8310cbbd2e1SMark F. Adams       /* get block */
8322e68589bSMark F. Adams       ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr);
8332e68589bSMark F. Adams       ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr);
8342e68589bSMark F. Adams       ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr);
8352e68589bSMark F. Adams       ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr);
8362e68589bSMark F. Adams       ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr);
8372e68589bSMark F. Adams 
8382e68589bSMark F. Adams       aggID = 0;
839e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(agg_llists,lid,&pos); CHKERRQ(ierr);
840e78576d6SMark F. Adams       while(pos){
841e78576d6SMark F. Adams         PetscInt gid1;
842ffc955d6SMark F. Adams         ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr);
843e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos); CHKERRQ(ierr);
844e78576d6SMark F. Adams 
8450cbbd2e1SMark F. Adams         if ( gid1 >= my0 && gid1 < Iend ) flid = gid1 - my0;
8460cbbd2e1SMark F. Adams         else {
8470cbbd2e1SMark F. Adams           ierr = GAMGTableFind( &fgid_flid, gid1, &flid ); CHKERRQ(ierr);
8480cbbd2e1SMark F. Adams           assert(flid>=0);
8490cbbd2e1SMark F. Adams         }
8502e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
85165d7b583SSatish Balay         data = &data_in[flid*bs];
8522e68589bSMark F. Adams         for ( kk = ii = 0; ii < bs ; ii++ ) {
8532e68589bSMark F. Adams           for ( jj = 0; jj < N ; jj++ ) {
8540cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8550cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8562e68589bSMark F. Adams           }
8572e68589bSMark F. Adams         }
8589057884aSMark F. Adams #ifdef OUT_AGGS
859b2a4f308SMark F. Adams         if (llev==1) {
8609057884aSMark F. Adams           char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8610cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
8620cbbd2e1SMark F. Adams           str[12] = col[(clid+17*mype)%6]; str[13] = sim[(clid+17*mype)%11];
863f7620de1SMatthew G Knepley           MM = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8640cbbd2e1SMark F. Adams           pj = gid1/MM; pi = gid1%MM;
865b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
866b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8679057884aSMark F. Adams         }
8689057884aSMark F. Adams #endif
8692e68589bSMark F. Adams         /* set fine IDs */
8702e68589bSMark F. Adams         for (kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8712e68589bSMark F. Adams 
8722e68589bSMark F. Adams         aggID++;
8730cbbd2e1SMark F. Adams       }
8742e68589bSMark F. Adams 
8752e68589bSMark F. Adams       /* pad with zeros */
8762e68589bSMark F. Adams       for ( ii = asz*bs; ii < Mdata ; ii++ ) {
8772e68589bSMark F. Adams 	for ( jj = 0; jj < N ; jj++, kk++ ) {
8782e68589bSMark F. Adams 	  qqc[jj*Mdata + ii] = .0;
8792e68589bSMark F. Adams 	}
8802e68589bSMark F. Adams       }
8812e68589bSMark F. Adams 
8822e68589bSMark F. Adams       ndone += aggID;
8832e68589bSMark F. Adams       /* QR */
88484df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8852e68589bSMark F. Adams       LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
88684df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
887c666adbfSMark F. Adams       if ( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRS error");
8882e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8892e68589bSMark F. Adams       {
8902e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8912e68589bSMark F. Adams         for ( jj = 0; jj < nSAvec ; jj++ ) {
8922e68589bSMark F. Adams           for ( ii = 0; ii < nSAvec ; ii++ ) {
8930cbbd2e1SMark F. Adams             assert(data[jj*out_data_stride + ii] == 1.e300);
8940cbbd2e1SMark F. Adams             if ( ii <= jj ) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8950cbbd2e1SMark F. Adams 	    else data[jj*out_data_stride + ii] = 0.;
8962e68589bSMark F. Adams           }
8972e68589bSMark F. Adams         }
8982e68589bSMark F. Adams       }
8992e68589bSMark F. Adams 
9002e68589bSMark F. Adams       /* get Q - row oriented */
9012e68589bSMark F. Adams       LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
902c666adbfSMark F. Adams       if ( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
9032e68589bSMark F. Adams 
9042e68589bSMark F. Adams       for ( ii = 0 ; ii < M ; ii++ ){
9052e68589bSMark F. Adams         for ( jj = 0 ; jj < N ; jj++ ) {
9062e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
9072e68589bSMark F. Adams         }
9082e68589bSMark F. Adams       }
9092e68589bSMark F. Adams 
9102e68589bSMark F. Adams       /* add diagonal block of P0 */
911c8b0795cSMark F. Adams       for (kk=0;kk<N;kk++) {
912c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
913c8b0795cSMark F. Adams       }
9142e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr);
9152e68589bSMark F. Adams 
9162e68589bSMark F. Adams       ierr = PetscFree( qqc );  CHKERRQ(ierr);
9172e68589bSMark F. Adams       ierr = PetscFree( qqr );  CHKERRQ(ierr);
9182e68589bSMark F. Adams       ierr = PetscFree( TAU );  CHKERRQ(ierr);
9192e68589bSMark F. Adams       ierr = PetscFree( WORK );  CHKERRQ(ierr);
9202e68589bSMark F. Adams       ierr = PetscFree( fids );  CHKERRQ(ierr);
921b43b03e9SMark F. Adams       clid++;
9220cbbd2e1SMark F. Adams     } /* coarse agg */
9230cbbd2e1SMark F. Adams   } /* for all fine nodes */
9240cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9250cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9262e68589bSMark F. Adams 
927c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */
9282e68589bSMark F. Adams /* MatGetSize( a_Prol, &kk, &jj ); */
929c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */
930e78576d6SMark 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); */
9312e68589bSMark F. Adams 
9329057884aSMark F. Adams #ifdef OUT_AGGS
933b2a4f308SMark F. Adams   if (llev==1) fclose(file);
9349057884aSMark F. Adams #endif
9350cbbd2e1SMark F. Adams   ierr = GAMGTableDestroy( &fgid_flid ); CHKERRQ(ierr);
9362e68589bSMark F. Adams 
9372e68589bSMark F. Adams   PetscFunctionReturn(0);
9382e68589bSMark F. Adams }
9392e68589bSMark F. Adams 
9402e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
9412e68589bSMark F. Adams /*
942c8b0795cSMark F. Adams    PCGAMGgraph_AGG
9432e68589bSMark F. Adams 
9442e68589bSMark F. Adams   Input Parameter:
9452e68589bSMark F. Adams    . pc - this
9462e68589bSMark F. Adams    . Amat - matrix on this fine level
9472e68589bSMark F. Adams   Output Parameter:
948c8b0795cSMark F. Adams    . a_Gmat -
9492e68589bSMark F. Adams */
9502e68589bSMark F. Adams #undef __FUNCT__
951c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG"
952c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_AGG( PC pc,
9532e68589bSMark F. Adams                                 const Mat Amat,
954c8b0795cSMark F. Adams                                 Mat *a_Gmat
955c8b0795cSMark F. Adams                                 )
956c8b0795cSMark F. Adams {
957c8b0795cSMark F. Adams   PetscErrorCode ierr;
958c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
959c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
960c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
961c8b0795cSMark F. Adams   const PetscReal vfilter = pc_gamg->threshold;
962c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
963c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
964e0940f08SMark F. Adams   Mat            Gmat;
965c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
96667744fedSMark F. Adams   PetscBool  /* set,flg , */symm;
967c8b0795cSMark F. Adams 
968c8b0795cSMark F. Adams   PetscFunctionBegin;
9690cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9700cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9710cbbd2e1SMark F. Adams #endif
972c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
973c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
974c8b0795cSMark F. Adams 
97567744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg); CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
976c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9770cbbd2e1SMark F. Adams 
9782d7fac45SMark F. Adams   ierr  = PCGAMGCreateGraph( Amat, &Gmat ); CHKERRQ( ierr );
9792d7fac45SMark F. Adams   ierr  = PCGAMGFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr );
980c8b0795cSMark F. Adams 
981e0940f08SMark F. Adams   *a_Gmat = Gmat;
982c8b0795cSMark F. Adams 
9830cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
9840cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
9850cbbd2e1SMark F. Adams #endif
986c8b0795cSMark F. Adams   PetscFunctionReturn(0);
987c8b0795cSMark F. Adams }
988c8b0795cSMark F. Adams 
989c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
990c8b0795cSMark F. Adams /*
991b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
992c8b0795cSMark F. Adams 
993c8b0795cSMark F. Adams   Input Parameter:
994e0940f08SMark F. Adams    . a_pc - this
995e0940f08SMark F. Adams   Input/Output Parameter:
9960cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
997c8b0795cSMark F. Adams   Output Parameter:
9980cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
999c8b0795cSMark F. Adams */
1000c8b0795cSMark F. Adams #undef __FUNCT__
1001b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
1002e0940f08SMark F. Adams PetscErrorCode PCGAMGCoarsen_AGG( PC a_pc,
1003e0940f08SMark F. Adams                                   Mat *a_Gmat1,
10040cbbd2e1SMark F. Adams                                   PetscCoarsenData **agg_lists
1005c8b0795cSMark F. Adams                                   )
1006c8b0795cSMark F. Adams {
1007c8b0795cSMark F. Adams   PetscErrorCode ierr;
1008e0940f08SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
1009c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1010c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
10110cbbd2e1SMark F. Adams   Mat             mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */
10129f3f12c8SMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
10130cbbd2e1SMark F. Adams   IS              perm;
1014c8b0795cSMark F. Adams   PetscInt        Ii,nloc,bs,n,m;
1015c8b0795cSMark F. Adams   PetscInt *permute;
1016c8b0795cSMark F. Adams   PetscBool *bIndexSet;
1017b43b03e9SMark F. Adams   MatCoarsen crs;
1018e0940f08SMark F. Adams   MPI_Comm        wcomm = ((PetscObject)Gmat1)->comm;
10190cbbd2e1SMark F. Adams   PetscMPIInt     mype,npe;
1020c8b0795cSMark F. Adams 
1021c8b0795cSMark F. Adams   PetscFunctionBegin;
10220cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
10230cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
10240cbbd2e1SMark F. Adams #endif
10250cbbd2e1SMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
10260cbbd2e1SMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
1027e0940f08SMark F. Adams   ierr = MatGetLocalSize( Gmat1, &n, &m ); CHKERRQ(ierr);
1028e0940f08SMark F. Adams   ierr = MatGetBlockSize( Gmat1, &bs ); CHKERRQ(ierr); assert(bs==1);
1029c8b0795cSMark F. Adams   nloc = n/bs;
1030c8b0795cSMark F. Adams 
1031e0940f08SMark F. Adams   if ( pc_gamg_agg->square_graph ) {
10329f3f12c8SMark F. Adams     if ( verbose > 1 ) {
10339f3f12c8SMark F. Adams       PetscPrintf(wcomm,"[%d]%s square graph\n",mype,__FUNCT__);
10349f3f12c8SMark F. Adams     }
1035878e152fSMark F. Adams     /* ierr = MatMatTransposeMult( Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); */
1036e0940f08SMark F. Adams     ierr = MatTransposeMatMult( Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 );
1037e0940f08SMark F. Adams     CHKERRQ(ierr);
10389f3f12c8SMark F. Adams     if ( verbose > 2 ) {
10399f3f12c8SMark F. Adams       PetscPrintf(wcomm,"[%d]%s square graph done\n",mype,__FUNCT__);
10409f3f12c8SMark F. Adams     }
1041e0940f08SMark F. Adams   }
1042e0940f08SMark F. Adams   else Gmat2 = Gmat1;
1043c8b0795cSMark F. Adams 
1044c8b0795cSMark F. Adams   /* get MIS aggs */
1045c8b0795cSMark F. Adams   /* randomize */
1046c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
1047c8b0795cSMark F. Adams   ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
1048c8b0795cSMark F. Adams   for ( Ii = 0; Ii < nloc ; Ii++ ){
1049c8b0795cSMark F. Adams     bIndexSet[Ii] = PETSC_FALSE;
1050c8b0795cSMark F. Adams     permute[Ii] = Ii;
1051c8b0795cSMark F. Adams   }
1052c8b0795cSMark F. Adams   srand(1); /* make deterministic */
1053c8b0795cSMark F. Adams   for ( Ii = 0; Ii < nloc ; Ii++ ) {
1054c8b0795cSMark F. Adams     PetscInt iSwapIndex = rand()%nloc;
1055c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1056c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1057c8b0795cSMark F. Adams       permute[iSwapIndex] = permute[Ii];
1058c8b0795cSMark F. Adams       permute[Ii] = iTemp;
1059c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1060c8b0795cSMark F. Adams     }
1061c8b0795cSMark F. Adams   }
1062c8b0795cSMark F. Adams   ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
1063c8b0795cSMark F. Adams 
10649f3f12c8SMark F. Adams   if ( verbose > 1 ) {
10659f3f12c8SMark F. Adams     PetscPrintf(wcomm,"[%d]%s coarsen graph\n",mype,__FUNCT__);
10669f3f12c8SMark F. Adams   }
10679f3f12c8SMark F. Adams 
1068c8b0795cSMark F. Adams   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);
1069c8b0795cSMark F. Adams   CHKERRQ(ierr);
10700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10710cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1072b43b03e9SMark F. Adams #endif
1073b43b03e9SMark F. Adams   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
10749057884aSMark F. Adams   /* ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); */
10759057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions( crs ); CHKERRQ(ierr);
1076b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
1077b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr);
1078b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
1079b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr);
1080b43b03e9SMark F. Adams   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
10810cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData( crs, agg_lists ); CHKERRQ(ierr); /* output */
1082b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
1083b43b03e9SMark F. Adams 
1084c8b0795cSMark F. Adams   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
1085c8b0795cSMark F. Adams   ierr = PetscFree( permute );  CHKERRQ(ierr);
10860cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10870cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1088b43b03e9SMark F. Adams #endif
10899f3f12c8SMark F. Adams   if ( verbose > 2 ) {
10909f3f12c8SMark F. Adams     PetscPrintf(wcomm,"[%d]%s coarsen graph done\n",mype,__FUNCT__);
10919f3f12c8SMark F. Adams   }
10929f3f12c8SMark F. Adams 
1093c8b0795cSMark F. Adams   /* smooth aggs */
1094e0940f08SMark F. Adams   if ( Gmat2 != Gmat1 ) {
10950cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10960cbbd2e1SMark F. Adams     ierr = smoothAggs( Gmat2, Gmat1, *agg_lists ); CHKERRQ(ierr);
1097c8b0795cSMark F. Adams     ierr = MatDestroy( &Gmat1 );  CHKERRQ(ierr);
1098e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
109941b27cdeSMark F. Adams     ierr = PetscCDGetMat( llist, &mat );  CHKERRQ(ierr);
11000cbbd2e1SMark F. Adams     if (mat) SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1101ef4ad70eSMark F. Adams   }
11020cbbd2e1SMark F. Adams   else {
11030cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
11040cbbd2e1SMark F. Adams     /* see if we have a matrix that takes pecedence (returned from MatCoarsenAppply) */
110541b27cdeSMark F. Adams     ierr = PetscCDGetMat( llist, &mat );   CHKERRQ(ierr);
11060cbbd2e1SMark F. Adams     if ( mat ) {
11070cbbd2e1SMark F. Adams       ierr = MatDestroy( &Gmat1 );  CHKERRQ(ierr);
11080cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
11090cbbd2e1SMark F. Adams     }
11100cbbd2e1SMark F. Adams   }
11110cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11120cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
11130cbbd2e1SMark F. Adams #endif
11149f3f12c8SMark F. Adams   if ( verbose > 2 ) {
11159f3f12c8SMark F. Adams     PetscPrintf(wcomm,"[%d]%s PCGAMGCoarsen_AGG done\n",mype,__FUNCT__);
11169f3f12c8SMark F. Adams   }
11179f3f12c8SMark F. Adams 
1118c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1119c8b0795cSMark F. Adams }
1120c8b0795cSMark F. Adams 
1121c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1122c8b0795cSMark F. Adams /*
11230cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1124c8b0795cSMark F. Adams 
1125c8b0795cSMark F. Adams  Input Parameter:
1126c8b0795cSMark F. Adams  . pc - this
1127c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1128c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
11290cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1130c8b0795cSMark F. Adams  Output Parameter:
1131c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1132c8b0795cSMark F. Adams  */
1133c8b0795cSMark F. Adams #undef __FUNCT__
11340cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
11350cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_AGG( PC pc,
1136c8b0795cSMark F. Adams                                       const Mat Amat,
1137c8b0795cSMark F. Adams                                       const Mat Gmat,
11380cbbd2e1SMark F. Adams                                       PetscCoarsenData *agg_lists,
1139c8b0795cSMark F. Adams                                       Mat *a_P_out
11402e68589bSMark F. Adams                                       )
11412e68589bSMark F. Adams {
11422e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
11432e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
11442e68589bSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
1145c8b0795cSMark F. Adams   const PetscInt data_cols = pc_gamg->data_cell_cols;
11462e68589bSMark F. Adams   PetscErrorCode ierr;
1147c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1148c8b0795cSMark F. Adams   Mat            Prol;
11492e68589bSMark F. Adams   PetscMPIInt    mype, npe;
11502e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
11510cbbd2e1SMark F. Adams   const PetscInt col_bs = data_cols;
1152c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1153c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
11542e68589bSMark F. Adams 
11552e68589bSMark F. Adams   PetscFunctionBegin;
11560cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
11570cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11580cbbd2e1SMark F. Adams #endif
11592e68589bSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
11602e68589bSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
11612e68589bSMark F. Adams   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
1162c8b0795cSMark F. Adams   ierr  = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr );
1163c8b0795cSMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
11642e68589bSMark F. Adams 
11652e68589bSMark F. Adams   /* get 'nLocalSelected' */
11660cbbd2e1SMark F. Adams   for ( ii=0, nLocalSelected = 0 ; ii < nloc ; ii++ ){
1167e78576d6SMark F. Adams     PetscBool ise;
11680cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1169e78576d6SMark F. Adams     ierr = PetscCDEmptyAt( agg_lists, ii, &ise ); CHKERRQ(ierr);
1170e78576d6SMark F. Adams     if ( !ise ) nLocalSelected++;
11712e68589bSMark F. Adams   }
11722e68589bSMark F. Adams 
11732e68589bSMark F. Adams   /* create prolongator, create P matrix */
1174a2f3521dSMark F. Adams   ierr = MatCreate( wcomm, &Prol ); CHKERRQ(ierr);
1175a2f3521dSMark F. Adams   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);
11762e68589bSMark F. Adams   CHKERRQ(ierr);
1177a2f3521dSMark F. Adams   ierr = MatSetBlockSizes( Prol, bs, col_bs ); CHKERRQ(ierr);
1178a2f3521dSMark F. Adams   ierr = MatSetType( Prol, MATAIJ );   CHKERRQ(ierr);
1179a2f3521dSMark F. Adams   ierr = MatSeqAIJSetPreallocation( Prol, data_cols, PETSC_NULL);  CHKERRQ(ierr);
1180a2f3521dSMark F. Adams   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, PETSC_NULL,data_cols, PETSC_NULL);CHKERRQ(ierr);
1181a2f3521dSMark F. Adams   /* nloc*bs, nLocalSelected*col_bs, */
1182a2f3521dSMark F. Adams   /* PETSC_DETERMINE, PETSC_DETERMINE, */
1183a2f3521dSMark F. Adams   /* data_cols, PETSC_NULL, data_cols, PETSC_NULL, */
1184a2f3521dSMark F. Adams   /* &Prol ); */
11852e68589bSMark F. Adams 
11862e68589bSMark F. Adams   /* can get all points "removed" */
1187c8b0795cSMark F. Adams   ierr =  MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr);
1188c8b0795cSMark F. Adams   if ( ii==0 ) {
11899f3f12c8SMark F. Adams     if ( verbose > 0 ) {
1190c8b0795cSMark F. Adams       PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__);
11912e68589bSMark F. Adams     }
11922e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
11932e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
11942e68589bSMark F. Adams     PetscFunctionReturn(0);
11952e68589bSMark F. Adams   }
11969f3f12c8SMark F. Adams   if ( verbose > 0 ) {
1197e78576d6SMark F. Adams     PetscPrintf(wcomm,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/col_bs);
1198c8b0795cSMark F. Adams   }
1199c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr);
12000cbbd2e1SMark F. Adams 
12010cbbd2e1SMark F. Adams   assert((kk-myCrs0)%col_bs==0);
1202c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
12030cbbd2e1SMark F. Adams   assert((kk/col_bs-myCrs0)==nLocalSelected);
12042e68589bSMark F. Adams 
12052e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
12060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12070cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
12082e68589bSMark F. Adams #endif
1209c8b0795cSMark F. Adams   if (npe > 1) { /*  */
12102e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
12112e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr);
12122e68589bSMark F. Adams     for ( jj = 0 ; jj < data_cols ; jj++ ){
1213c8b0795cSMark F. Adams       for ( kk = 0 ; kk < bs ; kk++) {
1214a2f3521dSMark F. Adams         PetscInt ii,stride;
1215c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1216c8b0795cSMark F. Adams         for ( ii = 0 ; ii < nloc ; ii++, tp += bs ){
12172e68589bSMark F. Adams           tmp_ldata[ii] = *tp;
12182e68589bSMark F. Adams         }
1219a2f3521dSMark F. Adams         ierr = PCGAMGGetDataWithGhosts( Gmat, 1, tmp_ldata, &stride, &tmp_gdata );
12202e68589bSMark F. Adams         CHKERRQ(ierr);
1221a2f3521dSMark F. Adams 
12222e68589bSMark F. Adams         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1223a2f3521dSMark F. Adams           ierr = PetscMalloc( stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr);
1224a2f3521dSMark F. Adams           nbnodes = bs*stride;
12252e68589bSMark F. Adams         }
1226a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
1227a2f3521dSMark F. Adams         for ( ii = 0 ; ii < stride ; ii++, tp2 += bs ){
12282e68589bSMark F. Adams           *tp2 = tmp_gdata[ii];
12292e68589bSMark F. Adams         }
12302e68589bSMark F. Adams         ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr);
12312e68589bSMark F. Adams       }
12322e68589bSMark F. Adams     }
12332e68589bSMark F. Adams     ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr);
12342e68589bSMark F. Adams   }
12352e68589bSMark F. Adams   else {
1236c8b0795cSMark F. Adams     nbnodes = bs*nloc;
1237c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
12382e68589bSMark F. Adams   }
12392e68589bSMark F. Adams 
12402e68589bSMark F. Adams   /* get P0 */
12412e68589bSMark F. Adams   if ( npe > 1 ){
12422e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1243a2f3521dSMark F. Adams     PetscInt stride;
12442e68589bSMark F. Adams 
12452e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr);
12462e68589bSMark F. Adams     for (kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1247a2f3521dSMark F. Adams     ierr = PCGAMGGetDataWithGhosts( Gmat, 1, fid_glid_loc, &stride, &fiddata );
12482e68589bSMark F. Adams     CHKERRQ(ierr);
1249a2f3521dSMark F. Adams     ierr = PetscMalloc( stride*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
1250a2f3521dSMark F. Adams     for (kk=0;kk<stride;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
12512e68589bSMark F. Adams     ierr = PetscFree( fiddata ); CHKERRQ(ierr);
1252a2f3521dSMark F. Adams 
1253a2f3521dSMark F. Adams     assert(stride==nbnodes/bs);
12542e68589bSMark F. Adams     ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr);
12552e68589bSMark F. Adams   }
12562e68589bSMark F. Adams   else {
12572e68589bSMark F. Adams     ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
12582e68589bSMark F. Adams     for (kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk;
12592e68589bSMark F. Adams   }
12600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12610cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
12620cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
12632e68589bSMark F. Adams #endif
1264c8b0795cSMark F. Adams   {
1265ffc955d6SMark F. Adams     PetscReal *data_out = PETSC_NULL;
12660cbbd2e1SMark F. Adams     ierr = formProl0( agg_lists, bs, data_cols, myCrs0, nbnodes,
1267c8b0795cSMark F. Adams                       data_w_ghost, flid_fgid, &data_out, Prol );
12682e68589bSMark F. Adams     CHKERRQ(ierr);
1269c8b0795cSMark F. Adams     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
1270c8b0795cSMark F. Adams     pc_gamg->data = data_out;
1271c8b0795cSMark F. Adams     pc_gamg->data_cell_rows = data_cols;
1272c8b0795cSMark F. Adams     pc_gamg->data_sz = data_cols*data_cols*nLocalSelected;
1273c8b0795cSMark F. Adams   }
12740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12750cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1276c8b0795cSMark F. Adams #endif
12772e68589bSMark F. Adams   if (npe > 1) ierr = PetscFree( data_w_ghost );      CHKERRQ(ierr);
12782e68589bSMark F. Adams   ierr = PetscFree( flid_fgid ); CHKERRQ(ierr);
12792e68589bSMark F. Adams 
1280c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
12810cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
12820cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
12830cbbd2e1SMark F. Adams #endif
1284c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1285c8b0795cSMark F. Adams }
1286c8b0795cSMark F. Adams 
1287c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1288c8b0795cSMark F. Adams /*
12890cbbd2e1SMark F. Adams    PCGAMGOptprol_AGG
1290c8b0795cSMark F. Adams 
1291c8b0795cSMark F. Adams   Input Parameter:
1292c8b0795cSMark F. Adams    . pc - this
1293c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1294c8b0795cSMark F. Adams  In/Output Parameter:
1295c8b0795cSMark F. Adams    . a_P_out - prolongation operator to the next level
1296c8b0795cSMark F. Adams */
1297c8b0795cSMark F. Adams #undef __FUNCT__
12980cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG"
12990cbbd2e1SMark F. Adams PetscErrorCode PCGAMGOptprol_AGG( PC pc,
1300c8b0795cSMark F. Adams                                   const Mat Amat,
1301c8b0795cSMark F. Adams                                   Mat *a_P
1302c8b0795cSMark F. Adams                                   )
1303c8b0795cSMark F. Adams {
1304c8b0795cSMark F. Adams   PetscErrorCode ierr;
1305c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
1306c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1307c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
1308c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1309c8b0795cSMark F. Adams   PetscInt       jj;
1310c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
1311c8b0795cSMark F. Adams   Mat            Prol = *a_P;
1312c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1313c8b0795cSMark F. Adams 
1314c8b0795cSMark F. Adams   PetscFunctionBegin;
13150cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13160cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
13170cbbd2e1SMark F. Adams #endif
1318c8b0795cSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
1319c8b0795cSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
1320c8b0795cSMark F. Adams 
13212e68589bSMark F. Adams   /* smooth P0 */
1322c8b0795cSMark F. Adams   for ( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){
13232e68589bSMark F. Adams     Mat tMat;
13242e68589bSMark F. Adams     Vec diag;
13252e68589bSMark F. Adams     PetscReal alpha, emax, emin;
13260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
13270cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
13282e68589bSMark F. Adams #endif
13292e68589bSMark F. Adams     if ( jj == 0 ) {
13302e68589bSMark F. Adams       KSP eksp;
13312e68589bSMark F. Adams       Vec bb, xx;
13322e68589bSMark F. Adams       PC pc;
13332e68589bSMark F. Adams       ierr = MatGetVecs( Amat, &bb, 0 );         CHKERRQ(ierr);
13342e68589bSMark F. Adams       ierr = MatGetVecs( Amat, &xx, 0 );         CHKERRQ(ierr);
13352e68589bSMark F. Adams       {
13362e68589bSMark F. Adams         PetscRandom    rctx;
13372e68589bSMark F. Adams         ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
13382e68589bSMark F. Adams         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
13392e68589bSMark F. Adams         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
13402e68589bSMark F. Adams         ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
13412e68589bSMark F. Adams       }
1342*e696ed0bSMark F. Adams 
1343*e696ed0bSMark F. Adams       /* zeroing out BC rows -- needed for crazy matrices */
1344*e696ed0bSMark F. Adams       {
1345*e696ed0bSMark F. Adams         PetscInt Istart,Iend,ncols,jj,Ii;
1346*e696ed0bSMark F. Adams         PetscScalar zero = 0.0;
1347*e696ed0bSMark F. Adams         ierr = MatGetOwnershipRange( Amat, &Istart, &Iend );    CHKERRQ(ierr);
1348*e696ed0bSMark F. Adams         for ( Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++ ) {
1349*e696ed0bSMark F. Adams           ierr = MatGetRow(Amat,Ii,&ncols,0,0); CHKERRQ(ierr);
1350*e696ed0bSMark F. Adams           if( ncols <= 1 ) {
1351*e696ed0bSMark F. Adams             ierr = VecSetValues( bb, 1, &Ii, &zero, INSERT_VALUES );  CHKERRQ(ierr);
1352*e696ed0bSMark F. Adams           }
1353*e696ed0bSMark F. Adams           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0); CHKERRQ(ierr);
1354*e696ed0bSMark F. Adams         }
1355*e696ed0bSMark F. Adams         ierr = VecAssemblyBegin(bb); CHKERRQ(ierr);
1356*e696ed0bSMark F. Adams         ierr = VecAssemblyEnd(bb); CHKERRQ(ierr);
1357*e696ed0bSMark F. Adams       }
1358*e696ed0bSMark F. Adams 
13592e68589bSMark F. Adams       ierr = KSPCreate(wcomm,&eksp);                            CHKERRQ(ierr);
13602e68589bSMark F. Adams       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);
13612e68589bSMark F. Adams       CHKERRQ(ierr);
13622e68589bSMark F. Adams       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
1363c2436cacSMark F. Adams       ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_");         CHKERRQ(ierr);
1364c2436cacSMark F. Adams       ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
1365c2436cacSMark F. Adams 
1366c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE );    CHKERRQ(ierr);
1367c2436cacSMark F. Adams       ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN );
1368c2436cacSMark F. Adams       CHKERRQ( ierr );
13692e68589bSMark F. Adams       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE );        CHKERRQ(ierr);
13702e68589bSMark F. Adams 
1371c2436cacSMark F. Adams       ierr = KSPGetPC( eksp, &pc );                              CHKERRQ( ierr );
1372c2436cacSMark F. Adams       ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr);  /* smoother in smoothed agg. */
1373c2436cacSMark F. Adams 
13742e68589bSMark F. Adams       /* solve - keep stuff out of logging */
13752e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
13762e68589bSMark F. Adams       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
13772e68589bSMark F. Adams       ierr = KSPSolve( eksp, bb, xx );                              CHKERRQ(ierr);
13782e68589bSMark F. Adams       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
13792e68589bSMark F. Adams       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
13802e68589bSMark F. Adams 
13812e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
13829f3f12c8SMark F. Adams       if ( verbose > 0 ) {
1383c8b0795cSMark F. Adams         PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
13842e68589bSMark F. Adams                     __FUNCT__,emax,emin,PCJACOBI);
13852e68589bSMark F. Adams       }
13862e68589bSMark F. Adams       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
13872e68589bSMark F. Adams       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
13882e68589bSMark F. Adams       ierr = KSPDestroy( &eksp );     CHKERRQ(ierr);
13892e68589bSMark F. Adams 
13902e68589bSMark F. Adams       if ( pc_gamg->emax_id == -1 ) {
13912e68589bSMark F. Adams         ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id );
13922e68589bSMark F. Adams         assert(pc_gamg->emax_id != -1 );
13932e68589bSMark F. Adams       }
13942e68589bSMark F. Adams       ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr);
13952e68589bSMark F. Adams     }
13962e68589bSMark F. Adams 
13972e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
13982e68589bSMark F. Adams     ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat );   CHKERRQ(ierr);
13992e68589bSMark F. Adams     ierr = MatGetVecs( Amat, &diag, 0 );    CHKERRQ(ierr);
14002e68589bSMark F. Adams     ierr = MatGetDiagonal( Amat, diag );    CHKERRQ(ierr); /* effectively PCJACOBI */
14012e68589bSMark F. Adams     ierr = VecReciprocal( diag );         CHKERRQ(ierr);
14022e68589bSMark F. Adams     ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr);
14032e68589bSMark F. Adams     ierr = VecDestroy( &diag );           CHKERRQ(ierr);
1404b4ec6429SMark F. Adams     alpha = -1.4/emax;
14052e68589bSMark F. Adams     ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN );           CHKERRQ(ierr);
14062e68589bSMark F. Adams     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
14072e68589bSMark F. Adams     Prol = tMat;
14080cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14090cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
14102e68589bSMark F. Adams #endif
14112e68589bSMark F. Adams   }
14120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
14130cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
14140cbbd2e1SMark F. Adams #endif
1415c8b0795cSMark F. Adams   *a_P = Prol;
1416c8b0795cSMark F. Adams 
1417c8b0795cSMark F. Adams   PetscFunctionReturn(0);
14182e68589bSMark F. Adams }
14192e68589bSMark F. Adams 
1420c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1421c8b0795cSMark F. Adams /*
1422a2f3521dSMark F. Adams    PCGAMGKKTProl_AGG
1423a2f3521dSMark F. Adams 
1424a2f3521dSMark F. Adams   Input Parameter:
1425a2f3521dSMark F. Adams    . pc - this
1426a2f3521dSMark F. Adams    . Prol11 - matrix on this fine level
1427a2f3521dSMark F. Adams    . A21 - matrix on this fine level
1428a2f3521dSMark F. Adams  In/Output Parameter:
1429a2f3521dSMark F. Adams    . a_P22 - prolongation operator to the next level
1430a2f3521dSMark F. Adams */
1431a2f3521dSMark F. Adams #undef __FUNCT__
1432a2f3521dSMark F. Adams #define __FUNCT__ "PCGAMGKKTProl_AGG"
1433a2f3521dSMark F. Adams PetscErrorCode PCGAMGKKTProl_AGG( PC pc,
1434a2f3521dSMark F. Adams                                   const Mat Prol11,
1435a2f3521dSMark F. Adams                                   const Mat A21,
1436a2f3521dSMark F. Adams                                   Mat *a_P22
1437a2f3521dSMark F. Adams                                   )
1438a2f3521dSMark F. Adams {
1439a2f3521dSMark F. Adams   PetscErrorCode ierr;
1440a2f3521dSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
1441a2f3521dSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1442a2f3521dSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
1443a2f3521dSMark F. Adams   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1444a2f3521dSMark F. Adams   PetscMPIInt    mype,npe;
1445a2f3521dSMark F. Adams   Mat            Prol22,Tmat,Gmat;
1446a2f3521dSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)pc)->comm;
1447a2f3521dSMark F. Adams   PetscCoarsenData *agg_lists;
1448a2f3521dSMark F. Adams 
1449a2f3521dSMark F. Adams   PetscFunctionBegin;
1450a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1451a2f3521dSMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0); CHKERRQ(ierr);
1452a2f3521dSMark F. Adams #endif
1453a2f3521dSMark F. Adams   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
1454a2f3521dSMark F. Adams   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
1455a2f3521dSMark F. Adams 
1456a2f3521dSMark F. Adams   /* form C graph */
1457a2f3521dSMark F. Adams   ierr = MatMatMult( A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);   CHKERRQ(ierr);
1458a2f3521dSMark F. Adams   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat ); CHKERRQ(ierr);
1459a2f3521dSMark F. Adams   ierr = MatDestroy(&Tmat);      CHKERRQ(ierr);
1460a2f3521dSMark F. Adams   ierr  = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose); CHKERRQ(ierr);
1461a2f3521dSMark F. Adams 
1462a2f3521dSMark F. Adams   /* coarsen constraints */
1463a2f3521dSMark F. Adams   {
1464a2f3521dSMark F. Adams     MatCoarsen crs;
1465a2f3521dSMark F. Adams     ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
1466a2f3521dSMark F. Adams     ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr);
1467a2f3521dSMark F. Adams     ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr);
1468a2f3521dSMark F. Adams     ierr = MatCoarsenSetVerbose( crs, verbose ); CHKERRQ(ierr);
1469a2f3521dSMark F. Adams     ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr);
1470a2f3521dSMark F. Adams     ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
1471a2f3521dSMark F. Adams     ierr = MatCoarsenGetData( crs, &agg_lists ); CHKERRQ(ierr);
1472a2f3521dSMark F. Adams     ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
1473a2f3521dSMark F. Adams   }
1474a2f3521dSMark F. Adams 
1475a2f3521dSMark F. Adams   /* form simple prolongation 'Prol22' */
1476a2f3521dSMark F. Adams   {
1477a2f3521dSMark F. Adams     PetscInt ii,mm,clid,my0,nloc,nLocalSelected;
1478a2f3521dSMark F. Adams     PetscScalar val = 1.0;
1479a2f3521dSMark F. Adams     /* get 'nLocalSelected' */
1480a2f3521dSMark F. Adams     ierr = MatGetLocalSize( Gmat, &nloc, &ii ); CHKERRQ(ierr);
1481a2f3521dSMark F. Adams     for ( ii=0, nLocalSelected = 0 ; ii < nloc ; ii++ ){
1482a2f3521dSMark F. Adams       PetscBool ise;
1483a2f3521dSMark F. Adams       /* filter out singletons 0 or 1? */
1484a2f3521dSMark F. Adams       ierr = PetscCDEmptyAt( agg_lists, ii, &ise ); CHKERRQ(ierr);
1485a2f3521dSMark F. Adams       if ( !ise ) nLocalSelected++;
1486a2f3521dSMark F. Adams     }
1487a2f3521dSMark F. Adams 
1488a2f3521dSMark F. Adams     ierr = MatCreate(wcomm,&Prol22);CHKERRQ(ierr);
1489a2f3521dSMark F. Adams     ierr = MatSetSizes( Prol22,nloc, nLocalSelected,
1490a2f3521dSMark F. Adams                         PETSC_DETERMINE, PETSC_DETERMINE);
1491a2f3521dSMark F. Adams     CHKERRQ(ierr);
1492a2f3521dSMark F. Adams     ierr = MatSetType( Prol22, MATAIJ );       CHKERRQ(ierr);
1493a2f3521dSMark F. Adams     ierr = MatSeqAIJSetPreallocation(Prol22,1,PETSC_NULL); CHKERRQ(ierr);
1494a2f3521dSMark F. Adams     ierr = MatMPIAIJSetPreallocation(Prol22,1,PETSC_NULL,1,PETSC_NULL);
1495a2f3521dSMark F. Adams     CHKERRQ(ierr);
1496a2f3521dSMark F. Adams     /* ierr = MatCreateAIJ( wcomm, */
1497a2f3521dSMark F. Adams     /*                      nloc, nLocalSelected, */
1498a2f3521dSMark F. Adams     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
1499a2f3521dSMark F. Adams     /*                      1, PETSC_NULL, 1, PETSC_NULL, */
1500a2f3521dSMark F. Adams     /*                      &Prol22 ); */
1501a2f3521dSMark F. Adams 
1502a2f3521dSMark F. Adams     ierr = MatGetOwnershipRange( Prol22, &my0, &ii );    CHKERRQ(ierr);
1503a2f3521dSMark F. Adams     nloc = ii - my0;
1504a2f3521dSMark F. Adams 
1505a2f3521dSMark F. Adams     /* make aggregates */
1506a2f3521dSMark F. Adams     for ( mm = clid = 0 ; mm < nloc ; mm++ ){
1507a2f3521dSMark F. Adams       ierr = PetscCDSizeAt( agg_lists, mm, &ii ); CHKERRQ(ierr);
1508a2f3521dSMark F. Adams       if ( ii > 0 ) {
1509a2f3521dSMark F. Adams         PetscInt asz=ii,cgid=my0+clid,rids[1000];
1510a2f3521dSMark F. Adams         PetscCDPos pos;
1511c666adbfSMark F. Adams         if (asz>1000)SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1512a2f3521dSMark F. Adams         ii = 0;
1513a2f3521dSMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos); CHKERRQ(ierr);
1514a2f3521dSMark F. Adams         while(pos){
1515a2f3521dSMark F. Adams           PetscInt gid1;
1516a2f3521dSMark F. Adams           ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr);
1517a2f3521dSMark F. Adams           ierr = PetscCDGetNextPos(agg_lists,mm,&pos); CHKERRQ(ierr);
1518a2f3521dSMark F. Adams 
1519a2f3521dSMark F. Adams           rids[ii++] = gid1;
1520a2f3521dSMark F. Adams         }
1521a2f3521dSMark F. Adams         assert(ii==asz);
1522a2f3521dSMark F. Adams         /* add diagonal block of P0 */
1523a2f3521dSMark F. Adams         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES); CHKERRQ(ierr);
1524a2f3521dSMark F. Adams 
1525a2f3521dSMark F. Adams         clid++;
1526a2f3521dSMark F. Adams       } /* coarse agg */
1527a2f3521dSMark F. Adams     } /* for all fine nodes */
1528a2f3521dSMark F. Adams     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);  CHKERRQ(ierr);
1529a2f3521dSMark F. Adams     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);    CHKERRQ(ierr);
1530a2f3521dSMark F. Adams   }
1531a2f3521dSMark F. Adams 
1532a2f3521dSMark F. Adams   /* clean up */
1533a2f3521dSMark F. Adams   ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
1534a2f3521dSMark F. Adams   ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
1535a2f3521dSMark F. Adams #if defined PETSC_USE_LOG
1536a2f3521dSMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1537a2f3521dSMark F. Adams #endif
1538a2f3521dSMark F. Adams   *a_P22 = Prol22;
1539a2f3521dSMark F. Adams 
1540a2f3521dSMark F. Adams   PetscFunctionReturn(0);
1541a2f3521dSMark F. Adams }
1542a2f3521dSMark F. Adams 
1543a2f3521dSMark F. Adams /* -------------------------------------------------------------------------- */
1544a2f3521dSMark F. Adams /*
1545c8b0795cSMark F. Adams    PCCreateGAMG_AGG
15462e68589bSMark F. Adams 
1547c8b0795cSMark F. Adams   Input Parameter:
1548c8b0795cSMark F. Adams    . pc -
1549c8b0795cSMark F. Adams */
1550c8b0795cSMark F. Adams #undef __FUNCT__
1551c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1552c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG( PC pc )
1553c8b0795cSMark F. Adams {
1554c8b0795cSMark F. Adams   PetscErrorCode  ierr;
1555c8b0795cSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1556c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1557c8b0795cSMark F. Adams   PC_GAMG_AGG      *pc_gamg_agg;
15582e68589bSMark F. Adams 
1559c8b0795cSMark F. Adams   PetscFunctionBegin;
1560ef820a72SMark F. Adams 
1561ef820a72SMark F. Adams   if( pc_gamg->subctx ){
1562ef820a72SMark F. Adams     /* call base class */
1563ef820a72SMark F. Adams     ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr);
1564ef820a72SMark F. Adams   }
1565ef820a72SMark F. Adams 
1566c8b0795cSMark F. Adams   /* create sub context for SA */
1567c8b0795cSMark F. Adams   ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr);
1568c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1569c8b0795cSMark F. Adams 
1570c8b0795cSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1571c8b0795cSMark F. Adams   pc->ops->destroy        = PCDestroy_AGG;
1572c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1573c8b0795cSMark F. Adams 
1574c8b0795cSMark F. Adams   /* set internal function pointers */
1575c8b0795cSMark F. Adams   pc_gamg->graph = PCGAMGgraph_AGG;
1576b43b03e9SMark F. Adams   pc_gamg->coarsen = PCGAMGCoarsen_AGG;
15770cbbd2e1SMark F. Adams   pc_gamg->prolongator = PCGAMGProlongator_AGG;
15780cbbd2e1SMark F. Adams   pc_gamg->optprol = PCGAMGOptprol_AGG;
1579a2f3521dSMark F. Adams   pc_gamg->formkktprol = PCGAMGKKTProl_AGG;
1580c8b0795cSMark F. Adams 
1581c8b0795cSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_AGG;
1582c8b0795cSMark F. Adams 
1583c8b0795cSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1584c8b0795cSMark F. Adams                                             "PCSetCoordinates_C",
1585c8b0795cSMark F. Adams                                             "PCSetCoordinates_AGG",
1586c8b0795cSMark F. Adams                                             PCSetCoordinates_AGG);
15872e68589bSMark F. Adams   PetscFunctionReturn(0);
15882e68589bSMark F. Adams }
1589