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