xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 9c334d8fdb557fc53fd345d68cbb3545b09ccab8)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 
5 #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
6 #include <petsc/private/kspimpl.h>
7 
8 #include <petscblaslapack.h>
9 
10 typedef struct {
11   PetscInt  nsmooths;
12   PetscBool sym_graph;
13   PetscInt square_graph;
14 } PC_GAMG_AGG;
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PCGAMGSetNSmooths"
18 /*@
19    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
20 
21    Not Collective on PC
22 
23    Input Parameters:
24 .  pc - the preconditioner context
25 
26    Options Database Key:
27 .  -pc_gamg_agg_nsmooths
28 
29    Level: intermediate
30 
31    Concepts: Aggregation AMG preconditioner
32 
33 .seealso: ()
34 @*/
35 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
36 {
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
41   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "PCGAMGSetNSmooths_AGG"
47 static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
48 {
49   PC_MG       *mg          = (PC_MG*)pc->data;
50   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
51   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
52 
53   PetscFunctionBegin;
54   pc_gamg_agg->nsmooths = n;
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "PCGAMGSetSymGraph"
60 /*@
61    PCGAMGSetSymGraph -
62 
63    Not Collective on PC
64 
65    Input Parameters:
66 .  pc - the preconditioner context
67 
68    Options Database Key:
69 .  -pc_gamg_sym_graph
70 
71    Level: intermediate
72 
73    Concepts: Aggregation AMG preconditioner
74 
75 .seealso: ()
76 @*/
77 PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
78 {
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
83   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "PCGAMGSetSymGraph_AGG"
89 static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n)
90 {
91   PC_MG       *mg          = (PC_MG*)pc->data;
92   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
93   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
94 
95   PetscFunctionBegin;
96   pc_gamg_agg->sym_graph = n;
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "PCGAMGSetSquareGraph"
102 /*@
103    PCGAMGSetSquareGraph -
104 
105    Not Collective on PC
106 
107    Input Parameters:
108 .  pc - the preconditioner context
109 
110    Options Database Key:
111 .  -pc_gamg_square_graph
112 
113    Level: intermediate
114 
115    Concepts: Aggregation AMG preconditioner
116 
117 .seealso: ()
118 @*/
119 PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
120 {
121   PetscErrorCode ierr;
122 
123   PetscFunctionBegin;
124   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "PCGAMGSetSquareGraph_AGG"
131 static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
132 {
133   PC_MG       *mg          = (PC_MG*)pc->data;
134   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
135   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
136 
137   PetscFunctionBegin;
138   pc_gamg_agg->square_graph = n;
139   PetscFunctionReturn(0);
140 }
141 
142 /* -------------------------------------------------------------------------- */
143 /*
144    PCSetFromOptions_GAMG_AGG
145 
146   Input Parameter:
147    . pc -
148 */
149 #undef __FUNCT__
150 #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
151 static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
152 {
153   PetscErrorCode ierr;
154   PC_MG          *mg          = (PC_MG*)pc->data;
155   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
156   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
157 
158   PetscFunctionBegin;
159   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
160   {
161     /* -pc_gamg_agg_nsmooths */
162     pc_gamg_agg->nsmooths = 1;
163 
164     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL);CHKERRQ(ierr);
165 
166     /* -pc_gamg_sym_graph */
167     pc_gamg_agg->sym_graph = PETSC_FALSE;
168 
169     ierr = PetscOptionsBool("-pc_gamg_sym_graph","Set for asymmetric matrices","PCGAMGSetSymGraph",pc_gamg_agg->sym_graph,&pc_gamg_agg->sym_graph,NULL);CHKERRQ(ierr);
170 
171     /* -pc_gamg_square_graph */
172     pc_gamg_agg->square_graph = 1;
173 
174     ierr = PetscOptionsInt("-pc_gamg_square_graph","Number of levels to square graph for faster coarsening and lower coarse grid complexity","PCGAMGSetSquareGraph",pc_gamg_agg->square_graph,&pc_gamg_agg->square_graph,NULL);CHKERRQ(ierr);
175   }
176   ierr = PetscOptionsTail();CHKERRQ(ierr);
177   PetscFunctionReturn(0);
178 }
179 
180 /* -------------------------------------------------------------------------- */
181 /*
182    PCDestroy_AGG
183 
184   Input Parameter:
185    . pc -
186 */
187 #undef __FUNCT__
188 #define __FUNCT__ "PCDestroy_GAMG_AGG"
189 static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
190 {
191   PetscErrorCode ierr;
192   PC_MG          *mg          = (PC_MG*)pc->data;
193   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
194 
195   PetscFunctionBegin;
196   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
197   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 /* -------------------------------------------------------------------------- */
202 /*
203    PCSetCoordinates_AGG
204      - collective
205 
206    Input Parameter:
207    . pc - the preconditioner context
208    . ndm - dimesion of data (used for dof/vertex for Stokes)
209    . a_nloc - number of vertices local
210    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
211 */
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "PCSetCoordinates_AGG"
215 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
216 {
217   PC_MG          *mg      = (PC_MG*)pc->data;
218   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
219   PetscErrorCode ierr;
220   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
221   Mat            mat = pc->pmat;
222 
223   PetscFunctionBegin;
224   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
225   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
226   nloc = a_nloc;
227 
228   /* SA: null space vectors */
229   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
230   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
231   else if (coords) {
232     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
233     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
234     if (ndm != ndf) {
235       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);
236     }
237   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
238   pc_gamg->data_cell_rows = ndatarows = ndf;
239   if (pc_gamg->data_cell_cols <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols);
240   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
241 
242   /* create data - syntactic sugar that should be refactored at some point */
243   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
244     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
245     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
246   }
247   /* copy data in - column oriented */
248   for (kk=0; kk<nloc; kk++) {
249     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
250     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
251     if (pc_gamg->data_cell_cols==1) *data = 1.0;
252     else {
253       /* translational modes */
254       for (ii=0;ii<ndatarows;ii++) {
255         for (jj=0;jj<ndatarows;jj++) {
256           if (ii==jj)data[ii*M + jj] = 1.0;
257           else data[ii*M + jj] = 0.0;
258         }
259       }
260 
261       /* rotational modes */
262       if (coords) {
263         if (ndm == 2) {
264           data   += 2*M;
265           data[0] = -coords[2*kk+1];
266           data[1] =  coords[2*kk];
267         } else {
268           data   += 3*M;
269           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
270           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
271           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
272         }
273       }
274     }
275   }
276 
277   pc_gamg->data_sz = arrsz;
278   PetscFunctionReturn(0);
279 }
280 
281 typedef PetscInt NState;
282 static const NState NOT_DONE=-2;
283 static const NState DELETED =-1;
284 static const NState REMOVED =-3;
285 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
286 
287 /* -------------------------------------------------------------------------- */
288 /*
289    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
290      - AGG-MG specific: clears singletons out of 'selected_2'
291 
292    Input Parameter:
293    . Gmat_2 - glabal matrix of graph (data not defined)   base (squared) graph
294    . Gmat_1 - base graph to grab with                 base graph
295    Input/Output Parameter:
296    . aggs_2 - linked list of aggs with gids)
297 */
298 #undef __FUNCT__
299 #define __FUNCT__ "smoothAggs"
300 static PetscErrorCode smoothAggs(Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
301 {
302   PetscErrorCode ierr;
303   PetscBool      isMPI;
304   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
305   MPI_Comm       comm;
306   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
307   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
308   const PetscInt nloc      = Gmat_2->rmap->n;
309   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
310   PetscInt       *lid_cprowID_1;
311   NState         *lid_state;
312   Vec            ghost_par_orig2;
313 
314   PetscFunctionBegin;
315   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
316   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
317 
318   if (PETSC_FALSE) {
319     PetscViewer viewer; char fname[32]; static int llev=0;
320     sprintf(fname,"Gmat2_%d.m",llev++);
321     PetscViewerASCIIOpen(comm,fname,&viewer);
322     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
323     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
324     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
325     ierr = PetscViewerDestroy(&viewer);
326   }
327 
328   /* get submatrices */
329   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
330   if (isMPI) {
331     /* grab matrix objects */
332     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
333     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
334     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
335     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
336     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
337     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
338 
339     /* force compressed row storage for B matrix in AuxMat */
340     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
341 
342     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
343     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
344     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
345       PetscInt lid = matB_1->compressedrow.rindex[ix];
346       lid_cprowID_1[lid] = ix;
347     }
348   } else {
349     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
350     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
351     lid_cprowID_1 = NULL;
352   }
353   if (nloc>0) {
354     if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
355     if (!(!matB_1 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
356     if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
357     if (!(!matB_2 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
358   }
359   /* get state of locals and selected gid for deleted */
360   ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr);
361   for (lid = 0; lid < nloc; lid++) {
362     lid_parent_gid[lid] = -1.0;
363     lid_state[lid]      = DELETED;
364   }
365 
366   /* set lid_state */
367   for (lid = 0; lid < nloc; lid++) {
368     PetscCDPos pos;
369     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
370     if (pos) {
371       PetscInt gid1;
372 
373       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
374       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
375       lid_state[lid] = gid1;
376     }
377   }
378 
379   /* map local to selected local, DELETED means a ghost owns it */
380   for (lid=kk=0; lid<nloc; lid++) {
381     NState state = lid_state[lid];
382     if (IS_SELECTED(state)) {
383       PetscCDPos pos;
384       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
385       while (pos) {
386         PetscInt gid1;
387         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
388         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
389 
390         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
391       }
392     }
393   }
394   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
395   if (isMPI) {
396     Vec tempVec;
397     /* get 'cpcol_1_state' */
398     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
399     for (kk=0,j=my0; kk<nloc; kk++,j++) {
400       PetscScalar v = (PetscScalar)lid_state[kk];
401       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
402     }
403     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
404     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
405     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
406     ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
407     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
408     /* get 'cpcol_2_state' */
409     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
411     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
412     /* get 'cpcol_2_par_orig' */
413     for (kk=0,j=my0; kk<nloc; kk++,j++) {
414       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
415       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
416     }
417     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
418     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
419     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
420     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
421     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
422     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
423 
424     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
425   } /* ismpi */
426 
427   /* doit */
428   for (lid=0; lid<nloc; lid++) {
429     NState state = lid_state[lid];
430     if (IS_SELECTED(state)) {
431       /* steal locals */
432       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
433       idx = matA_1->j + ii[lid];
434       for (j=0; j<n; j++) {
435         PetscInt lidj   = idx[j], sgid;
436         NState   statej = lid_state[lidj];
437         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
438           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
439           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
440             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
441             PetscCDPos pos,last=NULL;
442             /* looking for local from local so id_llist_2 works */
443             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
444             while (pos) {
445               PetscInt gid;
446               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
447               if (gid == gidj) {
448                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
449                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
450                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
451                 hav  = 1;
452                 break;
453               } else last = pos;
454 
455               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
456             }
457             if (hav!=1) {
458               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
459               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
460             }
461           } else {            /* I'm stealing this local, owned by a ghost */
462             if (sgid != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-pc_gamg_sym_graph true' to symetrize the graph or '-pc_gamg_threshold -1.0' if the matrix is structurally symmetric.");
463             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
464           }
465         }
466       } /* local neighbors */
467     } else if (state == DELETED && lid_cprowID_1) {
468       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
469       /* see if I have a selected ghost neighbor that will steal me */
470       if ((ix=lid_cprowID_1[lid]) != -1) {
471         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
472         idx = matB_1->j + ii[ix];
473         for (j=0; j<n; j++) {
474           PetscInt cpid   = idx[j];
475           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
476           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
477             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
478             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
479               PetscInt   hav=0,oldslidj=sgidold-my0;
480               PetscCDPos pos,last=NULL;
481               /* remove from 'oldslidj' list */
482               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
483               while (pos) {
484                 PetscInt gid;
485                 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
486                 if (lid+my0 == gid) {
487                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
488                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
489                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
490                   /* ghost (PetscScalar)statej will add this later */
491                   hav = 1;
492                   break;
493                 } else last = pos;
494 
495                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
496               }
497               if (hav!=1) {
498                 if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
499                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
500               }
501             } else {
502               /* ghosts remove this later */
503             }
504           }
505         }
506       }
507     } /* selected/deleted */
508   } /* node loop */
509 
510   if (isMPI) {
511     PetscScalar   *cpcol_2_parent,*cpcol_2_gid;
512     Vec           tempVec,ghostgids2,ghostparents2;
513     PetscInt      cpid,nghost_2;
514     GAMGHashTable gid_cpid;
515 
516     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
517     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
518 
519     /* get 'cpcol_2_parent' */
520     for (kk=0,j=my0; kk<nloc; kk++,j++) {
521       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
522     }
523     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
524     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
525     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
526     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
528     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
529 
530     /* get 'cpcol_2_gid' */
531     for (kk=0,j=my0; kk<nloc; kk++,j++) {
532       PetscScalar v = (PetscScalar)j;
533       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
534     }
535     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
536     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
537     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
538     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
541 
542     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
543 
544     /* look for deleted ghosts and add to table */
545     ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
546     for (cpid = 0; cpid < nghost_2; cpid++) {
547       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
548       if (state==DELETED) {
549         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
550         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
551         if (sgid_old == -1 && sgid_new != -1) {
552           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
553           ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
554         }
555       }
556     }
557 
558     /* look for deleted ghosts and see if they moved - remove it */
559     for (lid=0; lid<nloc; lid++) {
560       NState state = lid_state[lid];
561       if (IS_SELECTED(state)) {
562         PetscCDPos pos,last=NULL;
563         /* look for deleted ghosts and see if they moved */
564         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
565         while (pos) {
566           PetscInt gid;
567           ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
568 
569           if (gid < my0 || gid >= Iend) {
570             ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
571             if (cpid != -1) {
572               /* a moved ghost - */
573               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
574               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
575             } else last = pos;
576           } else last = pos;
577 
578           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
579         } /* loop over list of deleted */
580       } /* selected */
581     }
582     ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr);
583 
584     /* look at ghosts, see if they changed - and it */
585     for (cpid = 0; cpid < nghost_2; cpid++) {
586       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
587       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
588         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
589         PetscInt   slid_new=sgid_new-my0,hav=0;
590         PetscCDPos pos;
591         /* search for this gid to see if I have it */
592         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
593         while (pos) {
594           PetscInt gidj;
595           ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr);
596           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
597 
598           if (gidj == gid) { hav = 1; break; }
599         }
600         if (hav != 1) {
601           /* insert 'flidj' into head of llist */
602           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
603         }
604       }
605     }
606 
607     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
608     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
609     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
610     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
611     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
612     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
613     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
614     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
615   }
616 
617   ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr);
618   PetscFunctionReturn(0);
619 }
620 
621 /* -------------------------------------------------------------------------- */
622 /*
623    PCSetData_AGG - called if data is not set with PCSetCoordinates.
624       Looks in Mat for near null space.
625       Does not work for Stokes
626 
627   Input Parameter:
628    . pc -
629    . a_A - matrix to get (near) null space out of.
630 */
631 #undef __FUNCT__
632 #define __FUNCT__ "PCSetData_AGG"
633 static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
634 {
635   PetscErrorCode ierr;
636   PC_MG          *mg      = (PC_MG*)pc->data;
637   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
638   MatNullSpace   mnull;
639 
640   PetscFunctionBegin;
641   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
642   if (!mnull) {
643     PetscInt bs,NN,MM;
644     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
645     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
646     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
647     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
648   } else {
649     PetscReal *nullvec;
650     PetscBool has_const;
651     PetscInt  i,j,mlocal,nvec,bs;
652     const Vec *vecs; const PetscScalar *v;
653     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
654     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
655     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
656     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
657     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
658     for (i=0; i<nvec; i++) {
659       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
660       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
661       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
662     }
663     pc_gamg->data           = nullvec;
664     pc_gamg->data_cell_cols = (nvec+!!has_const);
665 
666     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
667 
668     pc_gamg->data_cell_rows = bs;
669   }
670   PetscFunctionReturn(0);
671 }
672 
673 /* -------------------------------------------------------------------------- */
674 /*
675  formProl0
676 
677    Input Parameter:
678    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
679    . bs - row block size
680    . nSAvec - column bs of new P
681    . my0crs - global index of start of locals
682    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
683    . data_in[data_stride*nSAvec] - local data on fine grid
684    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
685   Output Parameter:
686    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
687    . a_Prol - prolongation operator
688 */
689 #undef __FUNCT__
690 #define __FUNCT__ "formProl0"
691 static PetscErrorCode formProl0(PetscCoarsenData *agg_llists,PetscInt bs,PetscInt nSAvec,PetscInt my0crs,PetscInt data_stride,PetscReal data_in[],const PetscInt flid_fgid[],PetscReal **a_data_out,Mat a_Prol)
692 {
693   PetscErrorCode ierr;
694   PetscInt       Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
695   MPI_Comm       comm;
696   PetscMPIInt    rank;
697   PetscReal      *out_data;
698   PetscCDPos     pos;
699   GAMGHashTable  fgid_flid;
700 
701 /* #define OUT_AGGS */
702 #if defined(OUT_AGGS)
703   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
704 #endif
705 
706   PetscFunctionBegin;
707   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
708   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
709   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
710   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
711   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %d must be divisible by bs %D",Iend,Istart,bs);
712   Iend   /= bs;
713   nghosts = data_stride/bs - nloc;
714 
715   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
716   for (kk=0; kk<nghosts; kk++) {
717     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
718   }
719 
720 #if defined(OUT_AGGS)
721   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
722   if (llev==1) file = fopen(fname,"w");
723   MatGetSize(a_Prol, &pM, &jj);
724 #endif
725 
726   /* count selected -- same as number of cols of P */
727   for (nSelected=mm=0; mm<nloc; mm++) {
728     PetscBool ise;
729     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
730     if (!ise) nSelected++;
731   }
732   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
733   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
734   if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec);
735 
736   /* aloc space for coarse point data (output) */
737   out_data_stride = nSelected*nSAvec;
738 
739   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
740   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
741   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
742 
743   /* find points and set prolongation */
744   minsz = 100;
745   ndone = 0;
746   for (mm = clid = 0; mm < nloc; mm++) {
747     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
748     if (jj > 0) {
749       const PetscInt lid = mm, cgid = my0crs + clid;
750       PetscInt       cids[100]; /* max bs */
751       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
752       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
753       PetscScalar    *qqc,*qqr,*TAU,*WORK;
754       PetscInt       *fids;
755       PetscReal      *data;
756       /* count agg */
757       if (asz<minsz) minsz = asz;
758 
759       /* get block */
760       ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr);
761 
762       aggID = 0;
763       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
764       while (pos) {
765         PetscInt gid1;
766         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
767         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
768 
769         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
770         else {
771           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
772           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
773         }
774         /* copy in B_i matrix - column oriented */
775         data = &data_in[flid*bs];
776         for (ii = 0; ii < bs; ii++) {
777           for (jj = 0; jj < N; jj++) {
778             PetscReal d = data[jj*data_stride + ii];
779             qqc[jj*Mdata + aggID*bs + ii] = d;
780           }
781         }
782 #if defined(OUT_AGGS)
783         if (llev==1) {
784           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
785           PetscInt MM,pi,pj;
786           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
787           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
788           pj      = gid1/MM; pi = gid1%MM;
789           fprintf(file,str,(double)pi,(double)pj);
790           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
791         }
792 #endif
793         /* set fine IDs */
794         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
795 
796         aggID++;
797       }
798 
799       /* pad with zeros */
800       for (ii = asz*bs; ii < Mdata; ii++) {
801         for (jj = 0; jj < N; jj++, kk++) {
802           qqc[jj*Mdata + ii] = .0;
803         }
804       }
805 
806       ndone += aggID;
807       /* QR */
808       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
809       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
810       ierr = PetscFPTrapPop();CHKERRQ(ierr);
811       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
812       /* get R - column oriented - output B_{i+1} */
813       {
814         PetscReal *data = &out_data[clid*nSAvec];
815         for (jj = 0; jj < nSAvec; jj++) {
816           for (ii = 0; ii < nSAvec; ii++) {
817             if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",PETSC_MAX_REAL);
818            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
819            else data[jj*out_data_stride + ii] = 0.;
820           }
821         }
822       }
823 
824       /* get Q - row oriented */
825       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
826       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
827 
828       for (ii = 0; ii < M; ii++) {
829         for (jj = 0; jj < N; jj++) {
830           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
831         }
832       }
833 
834       /* add diagonal block of P0 */
835       for (kk=0; kk<N; kk++) {
836         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
837       }
838       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
839 
840       ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr);
841       clid++;
842     } /* coarse agg */
843   } /* for all fine nodes */
844   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
845   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
846 
847 /* ierr = MPIU_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
848 /* MatGetSize(a_Prol, &kk, &jj); */
849 /* ierr = MPIU_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
850 /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */
851 
852 #if defined(OUT_AGGS)
853   if (llev==1) fclose(file);
854 #endif
855   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "PCView_GAMG_AGG"
861 static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
862 {
863   PetscErrorCode ierr;
864   PC_MG          *mg      = (PC_MG*)pc->data;
865   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
866   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
867 
868   PetscFunctionBegin;
869   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
870   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
871   PetscFunctionReturn(0);
872 }
873 
874 /* -------------------------------------------------------------------------- */
875 /*
876    PCGAMGGraph_AGG
877 
878   Input Parameter:
879    . pc - this
880    . Amat - matrix on this fine level
881   Output Parameter:
882    . a_Gmat -
883 */
884 #undef __FUNCT__
885 #define __FUNCT__ "PCGAMGGraph_AGG"
886 static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
887 {
888   PetscErrorCode            ierr;
889   PC_MG                     *mg          = (PC_MG*)pc->data;
890   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
891   const PetscReal           vfilter      = pc_gamg->threshold;
892   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
893   Mat                       Gmat;
894   MPI_Comm                  comm;
895   PetscBool /* set,flg , */ symm;
896 
897   PetscFunctionBegin;
898   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
899   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
900 
901   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
902   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
903 
904   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
905   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
906 
907   *a_Gmat = Gmat;
908 
909   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
910   PetscFunctionReturn(0);
911 }
912 
913 /* -------------------------------------------------------------------------- */
914 /*
915    PCGAMGCoarsen_AGG
916 
917   Input Parameter:
918    . a_pc - this
919   Input/Output Parameter:
920    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
921   Output Parameter:
922    . agg_lists - list of aggregates
923 */
924 #undef __FUNCT__
925 #define __FUNCT__ "PCGAMGCoarsen_AGG"
926 static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
927 {
928   PetscErrorCode ierr;
929   PC_MG          *mg          = (PC_MG*)a_pc->data;
930   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
931   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
932   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
933   IS             perm;
934   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
935   PetscInt       *permute;
936   PetscBool      *bIndexSet;
937   MatCoarsen     crs;
938   MPI_Comm       comm;
939   PetscMPIInt    rank;
940   PetscReal      rr;
941   PetscInt       iSwapIndex;
942 
943   PetscFunctionBegin;
944   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
945   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
946   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
947   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
948   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
949   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
950   nloc = n/bs;
951 
952   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
953     ierr = PetscInfo2(a_pc,"Square Graph on level %d of %d to square\n",pc_gamg->current_level+1,pc_gamg_agg->square_graph);CHKERRQ(ierr);
954     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
955   } else Gmat2 = Gmat1;
956 
957   /* get MIS aggs - randomize */
958   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
959   ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
960   for (Ii = 0; Ii < nloc; Ii++) {
961     permute[Ii]   = Ii;
962   }
963   ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr);
964   for (Ii = 0; Ii < nloc; Ii++) {
965     ierr = PetscRandomGetValueReal(pc_gamg->random,&rr);CHKERRQ(ierr);
966     iSwapIndex = (PetscInt) (rr*nloc);
967     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
968       PetscInt iTemp = permute[iSwapIndex];
969       permute[iSwapIndex]   = permute[Ii];
970       permute[Ii]           = iTemp;
971       bIndexSet[iSwapIndex] = PETSC_TRUE;
972     }
973   }
974   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
975 
976   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
977 #if defined PETSC_GAMG_USE_LOG
978   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
979 #endif
980   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
981   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
982   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
983   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
984   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
985   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
986   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
987   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
988 
989   ierr = ISDestroy(&perm);CHKERRQ(ierr);
990   ierr = PetscFree(permute);CHKERRQ(ierr);
991 #if defined PETSC_GAMG_USE_LOG
992   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
993 #endif
994 
995   /* smooth aggs */
996   if (Gmat2 != Gmat1) {
997     const PetscCoarsenData *llist = *agg_lists;
998     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
999     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1000     *a_Gmat1 = Gmat2; /* output */
1001     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
1002     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1003   } else {
1004     const PetscCoarsenData *llist = *agg_lists;
1005     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
1006     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
1007     if (mat) {
1008       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1009       *a_Gmat1 = mat; /* output */
1010     }
1011   }
1012   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 /* -------------------------------------------------------------------------- */
1017 /*
1018  PCGAMGProlongator_AGG
1019 
1020  Input Parameter:
1021  . pc - this
1022  . Amat - matrix on this fine level
1023  . Graph - used to get ghost data for nodes in
1024  . agg_lists - list of aggregates
1025  Output Parameter:
1026  . a_P_out - prolongation operator to the next level
1027  */
1028 #undef __FUNCT__
1029 #define __FUNCT__ "PCGAMGProlongator_AGG"
1030 static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
1031 {
1032   PC_MG          *mg       = (PC_MG*)pc->data;
1033   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
1034   const PetscInt data_cols = pc_gamg->data_cell_cols;
1035   PetscErrorCode ierr;
1036   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1037   Mat            Prol;
1038   PetscMPIInt    rank, size;
1039   MPI_Comm       comm;
1040   const PetscInt col_bs = data_cols;
1041   PetscReal      *data_w_ghost;
1042   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1043   MatType        mtype;
1044 
1045   PetscFunctionBegin;
1046   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1047   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1048   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1049   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1050   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1051   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
1052   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
1053   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
1054 
1055   /* get 'nLocalSelected' */
1056   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1057     PetscBool ise;
1058     /* filter out singletons 0 or 1? */
1059     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1060     if (!ise) nLocalSelected++;
1061   }
1062 
1063   /* create prolongator, create P matrix */
1064   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
1065   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1066   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1067   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1068   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
1069   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
1070   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1071   /* nloc*bs, nLocalSelected*col_bs, */
1072   /* PETSC_DETERMINE, PETSC_DETERMINE, */
1073   /* data_cols, NULL, data_cols, NULL, */
1074   /* &Prol); */
1075 
1076   /* can get all points "removed" */
1077   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1078   if (!ii) {
1079     ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr);
1080     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
1081     *a_P_out = NULL;  /* out */
1082     PetscFunctionReturn(0);
1083   }
1084   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1085   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
1086 
1087   if ((kk-myCrs0) % col_bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs);
1088   myCrs0 = myCrs0/col_bs;
1089   if ((kk/col_bs-myCrs0) != nLocalSelected) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected);
1090 
1091   /* create global vector of data in 'data_w_ghost' */
1092 #if defined PETSC_GAMG_USE_LOG
1093   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1094 #endif
1095   if (size > 1) { /*  */
1096     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1097     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
1098     for (jj = 0; jj < data_cols; jj++) {
1099       for (kk = 0; kk < bs; kk++) {
1100         PetscInt        ii,stride;
1101         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1102         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
1103 
1104         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1105 
1106         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
1107           ierr    = PetscMalloc1(stride*bs*data_cols, &data_w_ghost);CHKERRQ(ierr);
1108           nbnodes = bs*stride;
1109         }
1110         tp2 = data_w_ghost + jj*bs*stride + kk;
1111         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
1112         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
1113       }
1114     }
1115     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1116   } else {
1117     nbnodes      = bs*nloc;
1118     data_w_ghost = (PetscReal*)pc_gamg->data;
1119   }
1120 
1121   /* get P0 */
1122   if (size > 1) {
1123     PetscReal *fid_glid_loc,*fiddata;
1124     PetscInt  stride;
1125 
1126     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
1127     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1128     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1129     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1130     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
1131     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1132 
1133     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
1134     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1135   } else {
1136     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
1137     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
1138   }
1139 #if defined PETSC_GAMG_USE_LOG
1140   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1141   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1142 #endif
1143   {
1144     PetscReal *data_out = NULL;
1145     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1146     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
1147 
1148     pc_gamg->data           = data_out;
1149     pc_gamg->data_cell_rows = data_cols;
1150     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1151   }
1152 #if defined PETSC_GAMG_USE_LOG
1153   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1154 #endif
1155   if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);}
1156   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
1157 
1158   *a_P_out = Prol;  /* out */
1159 
1160   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1161   PetscFunctionReturn(0);
1162 }
1163 
1164 /* -------------------------------------------------------------------------- */
1165 /*
1166    PCGAMGOptProlongator_AGG
1167 
1168   Input Parameter:
1169    . pc - this
1170    . Amat - matrix on this fine level
1171  In/Output Parameter:
1172    . a_P_out - prolongation operator to the next level
1173 */
1174 #undef __FUNCT__
1175 #define __FUNCT__ "PCGAMGOptProlongator_AGG"
1176 static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1177 {
1178   PetscErrorCode ierr;
1179   PC_MG          *mg          = (PC_MG*)pc->data;
1180   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1181   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1182   PetscInt       jj;
1183   Mat            Prol  = *a_P;
1184   MPI_Comm       comm;
1185 
1186   PetscFunctionBegin;
1187   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1188   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1189 
1190   /* smooth P0 */
1191   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1192     Mat       tMat;
1193     Vec       diag;
1194     PetscReal alpha, emax, emin;
1195 #if defined PETSC_GAMG_USE_LOG
1196     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1197 #endif
1198     if (!jj) {
1199       KSP         eksp;
1200       Vec         bb, xx;
1201       PC          epc;
1202 
1203       ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
1204       ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
1205 
1206       ierr = VecSetRandom(bb,pc_gamg->random);CHKERRQ(ierr);
1207 
1208       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1209       ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr);
1210       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
1211       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1212       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1213       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1214       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1215 
1216       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1217       ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
1218       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
1219 
1220       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1221       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1222 
1223       /* solve - keep stuff out of logging */
1224       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
1225       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
1226       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
1227       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
1228       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
1229 
1230       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1231       ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
1232       ierr = VecDestroy(&xx);CHKERRQ(ierr);
1233       ierr = VecDestroy(&bb);CHKERRQ(ierr);
1234       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
1235     }
1236 
1237     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
1238     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
1239     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
1240     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
1241     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
1242     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
1243     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1244     alpha = -1.4/emax;
1245     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
1246     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
1247     Prol  = tMat;
1248 #if defined PETSC_GAMG_USE_LOG
1249     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1250 #endif
1251   }
1252   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1253   *a_P = Prol;
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 /* -------------------------------------------------------------------------- */
1258 /*
1259    PCCreateGAMG_AGG
1260 
1261   Input Parameter:
1262    . pc -
1263 */
1264 #undef __FUNCT__
1265 #define __FUNCT__ "PCCreateGAMG_AGG"
1266 PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1267 {
1268   PetscErrorCode ierr;
1269   PC_MG          *mg      = (PC_MG*)pc->data;
1270   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1271   PC_GAMG_AGG    *pc_gamg_agg;
1272 
1273   PetscFunctionBegin;
1274   /* create sub context for SA */
1275   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1276   pc_gamg->subctx = pc_gamg_agg;
1277 
1278   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1279   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1280   /* reset does not do anything; setup not virtual */
1281 
1282   /* set internal function pointers */
1283   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
1284   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
1285   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1286   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
1287   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
1288   pc_gamg->ops->view              = PCView_GAMG_AGG;
1289 
1290   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr);
1291   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr);
1292   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr);
1293   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
1294   PetscFunctionReturn(0);
1295 }
1296