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