xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 186d4ecdc2f8f3d4a2a0951ddba16cdc5cd22fa1)
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 
254   PetscFunctionBegin;
255   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
256   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
257   nloc = a_nloc;
258 
259   /* SA: null space vectors */
260   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
261   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
262   else if (coords) {
263     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
264     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
265     if (ndm != ndf) {
266       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);
267     }
268   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
269   pc_gamg->data_cell_rows = ndatarows = ndf;
270   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);
271   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
272 
273   /* create data - syntactic sugar that should be refactored at some point */
274   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
275     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
276     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
277   }
278   /* copy data in - column oriented */
279   for (kk=0; kk<nloc; kk++) {
280     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
281     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
282     if (pc_gamg->data_cell_cols==1) *data = 1.0;
283     else {
284       /* translational modes */
285       for (ii=0;ii<ndatarows;ii++) {
286         for (jj=0;jj<ndatarows;jj++) {
287           if (ii==jj)data[ii*M + jj] = 1.0;
288           else data[ii*M + jj] = 0.0;
289         }
290       }
291 
292       /* rotational modes */
293       if (coords) {
294         if (ndm == 2) {
295           data   += 2*M;
296           data[0] = -coords[2*kk+1];
297           data[1] =  coords[2*kk];
298         } else {
299           data   += 3*M;
300           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
301           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
302           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
303         }
304       }
305     }
306   }
307 
308   pc_gamg->data_sz = arrsz;
309   PetscFunctionReturn(0);
310 }
311 EXTERN_C_END
312 
313 typedef PetscInt NState;
314 static const NState NOT_DONE=-2;
315 static const NState DELETED =-1;
316 static const NState REMOVED =-3;
317 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
318 
319 /* -------------------------------------------------------------------------- */
320 /*
321    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
322      - AGG-MG specific: clears singletons out of 'selected_2'
323 
324    Input Parameter:
325    . Gmat_2 - glabal matrix of graph (data not defined)
326    . Gmat_1 - base graph to grab with
327    Input/Output Parameter:
328    . aggs_2 - linked list of aggs with gids)
329 */
330 #undef __FUNCT__
331 #define __FUNCT__ "smoothAggs"
332 static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
333                                  const Mat Gmat_1,  /* base graph */
334                                  /* const IS selected_2, [nselected local] selected vertices */
335                                  PetscCoarsenData *aggs_2)  /* [nselected local] global ID of aggregate */
336 {
337   PetscErrorCode ierr;
338   PetscBool      isMPI;
339   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
340   MPI_Comm       comm;
341   PetscMPIInt    rank,size;
342   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
343   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
344   const PetscInt nloc      = Gmat_2->rmap->n;
345   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
346   PetscInt       *lid_cprowID_1;
347   NState         *lid_state;
348   Vec            ghost_par_orig2;
349 
350   PetscFunctionBegin;
351   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
352   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
353   ierr = MPI_Comm_size(comm, &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(comm,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       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 = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
755   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
756   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
757   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
758   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
759   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %d must be divisible by bs %D",Iend,Istart,bs);
760   Iend   /= bs;
761   nghosts = data_stride/bs - nloc;
762 
763   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
764   for (kk=0; kk<nghosts; kk++) {
765     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
766   }
767 
768 #if defined(OUT_AGGS)
769   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
770   if (llev==1) file = fopen(fname,"w");
771   MatGetSize(a_Prol, &pM, &jj);
772 #endif
773 
774   /* count selected -- same as number of cols of P */
775   for (nSelected=mm=0; mm<nloc; mm++) {
776     PetscBool ise;
777     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
778     if (!ise) nSelected++;
779   }
780   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
781   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
782   if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec);
783 
784   /* aloc space for coarse point data (output) */
785   out_data_stride = nSelected*nSAvec;
786 
787   ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr);
788   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=1.e300;
789   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
790 
791   /* find points and set prolongation */
792   minsz = 100;
793   ndone = 0;
794   for (mm = clid = 0; mm < nloc; mm++) {
795     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
796     if (jj > 0) {
797       const PetscInt lid = mm, cgid = my0crs + clid;
798       PetscInt       cids[100]; /* max bs */
799       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
800       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
801       PetscScalar    *qqc,*qqr,*TAU,*WORK;
802       PetscInt       *fids;
803       PetscReal      *data;
804       /* count agg */
805       if (asz<minsz) minsz = asz;
806 
807       /* get block */
808       ierr = PetscMalloc((Mdata*N)*sizeof(PetscScalar), &qqc);CHKERRQ(ierr);
809       ierr = PetscMalloc((M*N)*sizeof(PetscScalar), &qqr);CHKERRQ(ierr);
810       ierr = PetscMalloc(N*sizeof(PetscScalar), &TAU);CHKERRQ(ierr);
811       ierr = PetscMalloc(LWORK*sizeof(PetscScalar), &WORK);CHKERRQ(ierr);
812       ierr = PetscMalloc(M*sizeof(PetscInt), &fids);CHKERRQ(ierr);
813 
814       aggID = 0;
815       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
816       while (pos) {
817         PetscInt gid1;
818         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
819         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
820 
821         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
822         else {
823           ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
824           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
825         }
826         /* copy in B_i matrix - column oriented */
827         data = &data_in[flid*bs];
828         for (kk = ii = 0; ii < bs; ii++) {
829           for (jj = 0; jj < N; jj++) {
830             PetscReal d = data[jj*data_stride + ii];
831             qqc[jj*Mdata + aggID*bs + ii] = d;
832           }
833         }
834 #if defined(OUT_AGGS)
835         if (llev==1) {
836           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
837           PetscInt MM,pi,pj;
838           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
839           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
840           pj      = gid1/MM; pi = gid1%MM;
841           fprintf(file,str,(double)pi,(double)pj);
842           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
843         }
844 #endif
845         /* set fine IDs */
846         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
847 
848         aggID++;
849       }
850 
851       /* pad with zeros */
852       for (ii = asz*bs; ii < Mdata; ii++) {
853         for (jj = 0; jj < N; jj++, kk++) {
854           qqc[jj*Mdata + ii] = .0;
855         }
856       }
857 
858       ndone += aggID;
859       /* QR */
860       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
861       PetscStackCall("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
862       ierr = PetscFPTrapPop();CHKERRQ(ierr);
863       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
864       /* get R - column oriented - output B_{i+1} */
865       {
866         PetscReal *data = &out_data[clid*nSAvec];
867         for (jj = 0; jj < nSAvec; jj++) {
868           for (ii = 0; ii < nSAvec; ii++) {
869            if (data[jj*out_data_stride + ii] != 1.e300) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != 1.e300");
870            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
871            else data[jj*out_data_stride + ii] = 0.;
872           }
873         }
874       }
875 
876       /* get Q - row oriented */
877       PetscStackCall("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
878       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
879 
880       for (ii = 0; ii < M; ii++) {
881         for (jj = 0; jj < N; jj++) {
882           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
883         }
884       }
885 
886       /* add diagonal block of P0 */
887       for (kk=0; kk<N; kk++) {
888         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
889       }
890       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
891 
892       ierr = PetscFree(qqc);CHKERRQ(ierr);
893       ierr = PetscFree(qqr);CHKERRQ(ierr);
894       ierr = PetscFree(TAU);CHKERRQ(ierr);
895       ierr = PetscFree(WORK);CHKERRQ(ierr);
896       ierr = PetscFree(fids);CHKERRQ(ierr);
897       clid++;
898     } /* coarse agg */
899   } /* for all fine nodes */
900   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
901   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
902 
903 /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
904 /* MatGetSize(a_Prol, &kk, &jj); */
905 /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
906 /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */
907 
908 #if defined(OUT_AGGS)
909   if (llev==1) fclose(file);
910 #endif
911   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
912   PetscFunctionReturn(0);
913 }
914 
915 /* -------------------------------------------------------------------------- */
916 /*
917    PCGAMGgraph_AGG
918 
919   Input Parameter:
920    . pc - this
921    . Amat - matrix on this fine level
922   Output Parameter:
923    . a_Gmat -
924 */
925 #undef __FUNCT__
926 #define __FUNCT__ "PCGAMGgraph_AGG"
927 PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat)
928 {
929   PetscErrorCode            ierr;
930   PC_MG                     *mg          = (PC_MG*)pc->data;
931   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
932   const PetscInt            verbose      = pc_gamg->verbose;
933   const PetscReal           vfilter      = pc_gamg->threshold;
934   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
935   PetscMPIInt               rank,size;
936   Mat                       Gmat;
937   MPI_Comm                  comm;
938   PetscBool /* set,flg , */ symm;
939 
940   PetscFunctionBegin;
941   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
942 #if defined PETSC_USE_LOG
943   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
944 #endif
945   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
946   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
947 
948   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
949   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
950 
951   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
952   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
953 
954   *a_Gmat = Gmat;
955 
956 #if defined PETSC_USE_LOG
957   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
958 #endif
959   PetscFunctionReturn(0);
960 }
961 
962 /* -------------------------------------------------------------------------- */
963 /*
964    PCGAMGCoarsen_AGG
965 
966   Input Parameter:
967    . a_pc - this
968   Input/Output Parameter:
969    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
970   Output Parameter:
971    . agg_lists - list of aggregates
972 */
973 #undef __FUNCT__
974 #define __FUNCT__ "PCGAMGCoarsen_AGG"
975 PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
976 {
977   PetscErrorCode ierr;
978   PC_MG          *mg          = (PC_MG*)a_pc->data;
979   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
980   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
981   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
982   const PetscInt verbose = pc_gamg->verbose;
983   IS             perm;
984   PetscInt       Ii,nloc,bs,n,m;
985   PetscInt       *permute;
986   PetscBool      *bIndexSet;
987   MatCoarsen     crs;
988   MPI_Comm       comm;
989   PetscMPIInt    rank,size;
990 
991   PetscFunctionBegin;
992 #if defined PETSC_USE_LOG
993   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
994 #endif
995   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
996   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
997   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
998   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
999   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
1000   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
1001   nloc = n/bs;
1002 
1003   if (pc_gamg_agg->square_graph) {
1004     if (verbose > 1) PetscPrintf(comm,"[%d]%s square graph\n",rank,__FUNCT__);
1005     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
1006     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
1007     if (verbose > 2) {
1008       ierr = PetscPrintf(comm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr);
1009       /* check for symetry */
1010       if (verbose > 4) {
1011 
1012       }
1013     }
1014   } else Gmat2 = Gmat1;
1015 
1016   /* get MIS aggs */
1017   /* randomize */
1018   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
1019   ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
1020   for (Ii = 0; Ii < nloc; Ii++) {
1021     bIndexSet[Ii] = PETSC_FALSE;
1022     permute[Ii]   = Ii;
1023   }
1024   srand(1); /* make deterministic */
1025   for (Ii = 0; Ii < nloc; Ii++) {
1026     PetscInt iSwapIndex = rand()%nloc;
1027     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1028       PetscInt iTemp = permute[iSwapIndex];
1029       permute[iSwapIndex]   = permute[Ii];
1030       permute[Ii]           = iTemp;
1031       bIndexSet[iSwapIndex] = PETSC_TRUE;
1032     }
1033   }
1034   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
1035 
1036   if (verbose > 1) PetscPrintf(comm,"[%d]%s coarsen graph\n",rank,__FUNCT__);
1037 
1038   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
1039 #if defined PETSC_GAMG_USE_LOG
1040   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1041 #endif
1042   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
1043   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
1044   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1045   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1046   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1047   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
1048   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1049   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1050   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1051   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1052 
1053   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1054   ierr = PetscFree(permute);CHKERRQ(ierr);
1055 #if defined PETSC_GAMG_USE_LOG
1056   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1057 #endif
1058   if (verbose > 2) PetscPrintf(comm,"[%d]%s coarsen graph done\n",rank,__FUNCT__);
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(comm,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) PetscPrintf(comm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__);
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 /* -------------------------------------------------------------------------- */
1085 /*
1086  PCGAMGProlongator_AGG
1087 
1088  Input Parameter:
1089  . pc - this
1090  . Amat - matrix on this fine level
1091  . Graph - used to get ghost data for nodes in
1092  . agg_lists - list of aggregates
1093  Output Parameter:
1094  . a_P_out - prolongation operator to the next level
1095  */
1096 #undef __FUNCT__
1097 #define __FUNCT__ "PCGAMGProlongator_AGG"
1098 PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
1099 {
1100   PC_MG          *mg       = (PC_MG*)pc->data;
1101   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
1102   const PetscInt verbose   = pc_gamg->verbose;
1103   const PetscInt data_cols = pc_gamg->data_cell_cols;
1104   PetscErrorCode ierr;
1105   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1106   Mat            Prol;
1107   PetscMPIInt    rank, size;
1108   MPI_Comm       comm;
1109   const PetscInt col_bs = data_cols;
1110   PetscReal      *data_w_ghost;
1111   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1112 
1113   PetscFunctionBegin;
1114   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1115 #if defined PETSC_USE_LOG
1116   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1117 #endif
1118   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1119   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1120   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1121   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
1122   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
1123   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
1124 
1125   /* get 'nLocalSelected' */
1126   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1127     PetscBool ise;
1128     /* filter out singletons 0 or 1? */
1129     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1130     if (!ise) nLocalSelected++;
1131   }
1132 
1133   /* create prolongator, create P matrix */
1134   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1135   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1136   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1137   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
1138   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr);
1139   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr);
1140   /* nloc*bs, nLocalSelected*col_bs, */
1141   /* PETSC_DETERMINE, PETSC_DETERMINE, */
1142   /* data_cols, NULL, data_cols, NULL, */
1143   /* &Prol); */
1144 
1145   /* can get all points "removed" */
1146   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1147   if (ii==0) {
1148     if (verbose > 0) PetscPrintf(comm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__);
1149     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
1150     *a_P_out = NULL;  /* out */
1151     PetscFunctionReturn(0);
1152   }
1153   if (verbose > 0) PetscPrintf(comm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs);
1154   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
1155 
1156   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);
1157   myCrs0 = myCrs0/col_bs;
1158   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);
1159 
1160   /* create global vector of data in 'data_w_ghost' */
1161 #if defined PETSC_GAMG_USE_LOG
1162   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1163 #endif
1164   if (size > 1) { /*  */
1165     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1166     ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr);
1167     for (jj = 0; jj < data_cols; jj++) {
1168       for (kk = 0; kk < bs; kk++) {
1169         PetscInt        ii,stride;
1170         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1171         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
1172 
1173         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1174 
1175         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1176           ierr    = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr);
1177           nbnodes = bs*stride;
1178         }
1179         tp2 = data_w_ghost + jj*bs*stride + kk;
1180         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
1181         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
1182       }
1183     }
1184     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1185   } else {
1186     nbnodes      = bs*nloc;
1187     data_w_ghost = (PetscReal*)pc_gamg->data;
1188   }
1189 
1190   /* get P0 */
1191   if (size > 1) {
1192     PetscReal *fid_glid_loc,*fiddata;
1193     PetscInt  stride;
1194 
1195     ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr);
1196     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1197     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1198     ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1199     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
1200     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1201 
1202     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
1203     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1204   } else {
1205     ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1206     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
1207   }
1208 #if defined PETSC_GAMG_USE_LOG
1209   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1210   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1211 #endif
1212   {
1213     PetscReal *data_out = NULL;
1214     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1215     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
1216 
1217     pc_gamg->data           = data_out;
1218     pc_gamg->data_cell_rows = data_cols;
1219     pc_gamg->data_sz        = data_cols*data_cols*nLocalSelected;
1220   }
1221 #if defined PETSC_GAMG_USE_LOG
1222   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1223 #endif
1224   if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
1225   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
1226 
1227   *a_P_out = Prol;  /* out */
1228 #if defined PETSC_USE_LOG
1229   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1230 #endif
1231   PetscFunctionReturn(0);
1232 }
1233 
1234 /* -------------------------------------------------------------------------- */
1235 /*
1236    PCGAMGOptprol_AGG
1237 
1238   Input Parameter:
1239    . pc - this
1240    . Amat - matrix on this fine level
1241  In/Output Parameter:
1242    . a_P_out - prolongation operator to the next level
1243 */
1244 #undef __FUNCT__
1245 #define __FUNCT__ "PCGAMGOptprol_AGG"
1246 PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P)
1247 {
1248   PetscErrorCode ierr;
1249   PC_MG          *mg          = (PC_MG*)pc->data;
1250   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1251   const PetscInt verbose      = pc_gamg->verbose;
1252   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1253   PetscInt       jj;
1254   PetscMPIInt    rank,size;
1255   Mat            Prol  = *a_P;
1256   MPI_Comm       comm;
1257 
1258   PetscFunctionBegin;
1259   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1260 #if defined PETSC_USE_LOG
1261   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
1262 #endif
1263   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1264   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1265 
1266   /* smooth P0 */
1267   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1268     Mat       tMat;
1269     Vec       diag;
1270     PetscReal alpha, emax, emin;
1271 #if defined PETSC_GAMG_USE_LOG
1272     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1273 #endif
1274     if (jj == 0) {
1275       KSP eksp;
1276       Vec bb, xx;
1277       PC  epc;
1278       ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr);
1279       ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr);
1280       {
1281         PetscRandom rctx;
1282         ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
1283         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1284         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
1285         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1286       }
1287 
1288       /* zeroing out BC rows -- needed for crazy matrices */
1289       {
1290         PetscInt    Istart,Iend,ncols,jj,Ii;
1291         PetscScalar zero = 0.0;
1292         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1293         for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
1294           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1295           if (ncols <= 1) {
1296             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1297           }
1298           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1299         }
1300         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1301         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1302       }
1303 
1304       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1305       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
1306       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1307       ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1308       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1309       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1310 
1311       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1312       ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1313       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
1314 
1315       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1316       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1317 
1318       /* solve - keep stuff out of logging */
1319       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
1320       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
1321       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
1322       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
1323       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
1324 
1325       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1326       if (verbose > 0) {
1327         PetscPrintf(comm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
1328                     __FUNCT__,emax,emin,PCJACOBI);
1329       }
1330       ierr = VecDestroy(&xx);CHKERRQ(ierr);
1331       ierr = VecDestroy(&bb);CHKERRQ(ierr);
1332       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
1333 
1334       if (pc_gamg->emax_id == -1) {
1335         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
1336         if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1");
1337       }
1338       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
1339     }
1340 
1341     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
1342     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
1343     ierr  = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr);
1344     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
1345     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
1346     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
1347     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1348     alpha = -1.4/emax;
1349     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
1350     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
1351     Prol  = tMat;
1352 #if defined PETSC_GAMG_USE_LOG
1353     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1354 #endif
1355   }
1356 #if defined PETSC_USE_LOG
1357   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
1358 #endif
1359   *a_P = Prol;
1360   PetscFunctionReturn(0);
1361 }
1362 
1363 /* -------------------------------------------------------------------------- */
1364 /*
1365    PCGAMGKKTProl_AGG
1366 
1367   Input Parameter:
1368    . pc - this
1369    . Prol11 - matrix on this fine level
1370    . A21 - matrix on this fine level
1371  In/Output Parameter:
1372    . a_P22 - prolongation operator to the next level
1373 */
1374 #undef __FUNCT__
1375 #define __FUNCT__ "PCGAMGKKTProl_AGG"
1376 PetscErrorCode PCGAMGKKTProl_AGG(PC pc,const Mat Prol11,const Mat A21,Mat *a_P22)
1377 {
1378   PetscErrorCode   ierr;
1379   PC_MG            *mg      = (PC_MG*)pc->data;
1380   PC_GAMG          *pc_gamg = (PC_GAMG*)mg->innerctx;
1381   const PetscInt   verbose  = pc_gamg->verbose;
1382   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1383   PetscMPIInt      rank,size;
1384   Mat              Prol22,Tmat,Gmat;
1385   MPI_Comm         comm;
1386   PetscCoarsenData *agg_lists;
1387 
1388   PetscFunctionBegin;
1389   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1390 #if defined PETSC_USE_LOG
1391   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1392 #endif
1393   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1394   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1395 
1396   /* form C graph */
1397   ierr = MatMatMult(A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);CHKERRQ(ierr);
1398   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat);CHKERRQ(ierr);
1399   ierr = MatDestroy(&Tmat);CHKERRQ(ierr);
1400   ierr = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose);CHKERRQ(ierr);
1401 
1402   /* coarsen constraints */
1403   {
1404     MatCoarsen crs;
1405     ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
1406     ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
1407     ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
1408     ierr = MatCoarsenSetVerbose(crs, verbose);CHKERRQ(ierr);
1409     ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1410     ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1411     ierr = MatCoarsenGetData(crs, &agg_lists);CHKERRQ(ierr);
1412     ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1413   }
1414 
1415   /* form simple prolongation 'Prol22' */
1416   {
1417     PetscInt    ii,mm,clid,my0,nloc,nLocalSelected;
1418     PetscScalar val = 1.0;
1419     /* get 'nLocalSelected' */
1420     ierr = MatGetLocalSize(Gmat, &nloc, &ii);CHKERRQ(ierr);
1421     for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1422       PetscBool ise;
1423       /* filter out singletons 0 or 1? */
1424       ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1425       if (!ise) nLocalSelected++;
1426     }
1427 
1428     ierr = MatCreate(comm,&Prol22);CHKERRQ(ierr);
1429     ierr = MatSetSizes(Prol22,nloc, nLocalSelected,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1430     ierr = MatSetType(Prol22, MATAIJ);CHKERRQ(ierr);
1431     ierr = MatSeqAIJSetPreallocation(Prol22,1,NULL);CHKERRQ(ierr);
1432     ierr = MatMPIAIJSetPreallocation(Prol22,1,NULL,1,NULL);CHKERRQ(ierr);
1433     /* ierr = MatCreateAIJ(comm, */
1434     /*                      nloc, nLocalSelected, */
1435     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
1436     /*                      1, NULL, 1, NULL, */
1437     /*                      &Prol22); */
1438 
1439     ierr = MatGetOwnershipRange(Prol22, &my0, &ii);CHKERRQ(ierr);
1440     nloc = ii - my0;
1441 
1442     /* make aggregates */
1443     for (mm = clid = 0; mm < nloc; mm++) {
1444       ierr = PetscCDSizeAt(agg_lists, mm, &ii);CHKERRQ(ierr);
1445       if (ii > 0) {
1446         PetscInt   asz=ii,cgid=my0+clid,rids[1000];
1447         PetscCDPos pos;
1448         if (asz>1000) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1449         ii   = 0;
1450         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1451         while (pos) {
1452           PetscInt gid1;
1453           ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
1454           ierr = PetscCDGetNextPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1455 
1456           rids[ii++] = gid1;
1457         }
1458         if (ii != asz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D != asz %D",ii,asz);
1459         /* add diagonal block of P0 */
1460         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES);CHKERRQ(ierr);
1461 
1462         clid++;
1463       } /* coarse agg */
1464     } /* for all fine nodes */
1465     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1466     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1467   }
1468 
1469   /* clean up */
1470   ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
1471   ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
1472 #if defined PETSC_USE_LOG
1473   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1474 #endif
1475   *a_P22 = Prol22;
1476   PetscFunctionReturn(0);
1477 }
1478 
1479 /* -------------------------------------------------------------------------- */
1480 /*
1481    PCCreateGAMG_AGG
1482 
1483   Input Parameter:
1484    . pc -
1485 */
1486 #undef __FUNCT__
1487 #define __FUNCT__ "PCCreateGAMG_AGG"
1488 PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1489 {
1490   PetscErrorCode ierr;
1491   PC_MG          *mg      = (PC_MG*)pc->data;
1492   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1493   PC_GAMG_AGG    *pc_gamg_agg;
1494 
1495   PetscFunctionBegin;
1496   if (pc_gamg->subctx) {
1497     /* call base class */
1498     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
1499   }
1500 
1501   /* create sub context for SA */
1502   ierr            = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr);
1503   pc_gamg->subctx = pc_gamg_agg;
1504 
1505   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1506   pc->ops->destroy        = PCDestroy_AGG;
1507   /* reset does not do anything; setup not virtual */
1508 
1509   /* set internal function pointers */
1510   pc_gamg->graph       = PCGAMGgraph_AGG;
1511   pc_gamg->coarsen     = PCGAMGCoarsen_AGG;
1512   pc_gamg->prolongator = PCGAMGProlongator_AGG;
1513   pc_gamg->optprol     = PCGAMGOptprol_AGG;
1514   pc_gamg->formkktprol = PCGAMGKKTProl_AGG;
1515 
1516   pc_gamg->createdefaultdata = PCSetData_AGG;
1517 
1518   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_AGG",PCSetCoordinates_AGG);CHKERRQ(ierr);
1519   PetscFunctionReturn(0);
1520 }
1521