xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 3bf036e263fd78807e2931ff42d9ddcd8aae3fd4)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 
5 #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
6 #include <petsc-private/kspimpl.h>
7 
8 #include <assert.h>
9 #include <petscblaslapack.h>
10 
11 typedef struct {
12   PetscInt nsmooths;
13   PetscBool sym_graph;
14   PetscBool square_graph;
15 }PC_GAMG_AGG;
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "PCGAMGSetNSmooths"
19 /*@
20    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
21 
22    Not Collective on PC
23 
24    Input Parameters:
25 .  pc - the preconditioner context
26 
27    Options Database Key:
28 .  -pc_gamg_agg_nsmooths
29 
30    Level: intermediate
31 
32    Concepts: Aggregation AMG preconditioner
33 
34 .seealso: ()
35 @*/
36 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
37 {
38   PetscErrorCode ierr;
39 
40   PetscFunctionBegin;
41   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
42   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
43   PetscFunctionReturn(0);
44 }
45 
46 EXTERN_C_BEGIN
47 #undef __FUNCT__
48 #define __FUNCT__ "PCGAMGSetNSmooths_GAMG"
49 PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n)
50 {
51   PC_MG           *mg = (PC_MG*)pc->data;
52   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
53   PC_GAMG_AGG     *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
54 
55   PetscFunctionBegin;
56   pc_gamg_agg->nsmooths = n;
57   PetscFunctionReturn(0);
58 }
59 EXTERN_C_END
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "PCGAMGSetSymGraph"
63 /*@
64    PCGAMGSetSymGraph -
65 
66    Not Collective on PC
67 
68    Input Parameters:
69 .  pc - the preconditioner context
70 
71    Options Database Key:
72 .  -pc_gamg_sym_graph
73 
74    Level: intermediate
75 
76    Concepts: Aggregation AMG preconditioner
77 
78 .seealso: ()
79 @*/
80 PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
81 {
82   PetscErrorCode ierr;
83 
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 EXTERN_C_BEGIN
91 #undef __FUNCT__
92 #define __FUNCT__ "PCGAMGSetSymGraph_GAMG"
93 PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n)
94 {
95   PC_MG           *mg = (PC_MG*)pc->data;
96   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
97   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
98 
99   PetscFunctionBegin;
100   pc_gamg_agg->sym_graph = n;
101   PetscFunctionReturn(0);
102 }
103 EXTERN_C_END
104 
105 #undef __FUNCT__
106 #define __FUNCT__ "PCGAMGSetSquareGraph"
107 /*@
108    PCGAMGSetSquareGraph -
109 
110    Not Collective on PC
111 
112    Input Parameters:
113 .  pc - the preconditioner context
114 
115    Options Database Key:
116 .  -pc_gamg_square_graph
117 
118    Level: intermediate
119 
120    Concepts: Aggregation AMG preconditioner
121 
122 .seealso: ()
123 @*/
124 PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n)
125 {
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
130   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
131   PetscFunctionReturn(0);
132 }
133 
134 EXTERN_C_BEGIN
135 #undef __FUNCT__
136 #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG"
137 PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n)
138 {
139   PC_MG           *mg = (PC_MG*)pc->data;
140   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
141   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
142 
143   PetscFunctionBegin;
144   pc_gamg_agg->square_graph = n;
145   PetscFunctionReturn(0);
146 }
147 EXTERN_C_END
148 
149 /* -------------------------------------------------------------------------- */
150 /*
151    PCSetFromOptions_GAMG_AGG
152 
153   Input Parameter:
154    . pc -
155 */
156 #undef __FUNCT__
157 #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
158 PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc)
159 {
160   PetscErrorCode  ierr;
161   PC_MG           *mg = (PC_MG*)pc->data;
162   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
163   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
164   PetscBool        flag;
165 
166   PetscFunctionBegin;
167   /* call base class */
168   ierr = PCSetFromOptions_GAMG(pc);CHKERRQ(ierr);
169 
170   ierr = PetscOptionsHead("GAMG-AGG options");CHKERRQ(ierr);
171   {
172     /* -pc_gamg_agg_nsmooths */
173     pc_gamg_agg->nsmooths = 0;
174     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
175                            "smoothing steps for smoothed aggregation, usually 1 (0)",
176                            "PCGAMGSetNSmooths",
177                            pc_gamg_agg->nsmooths,
178                            &pc_gamg_agg->nsmooths,
179                            &flag);CHKERRQ(ierr);
180 
181     /* -pc_gamg_sym_graph */
182     pc_gamg_agg->sym_graph = PETSC_FALSE;
183     ierr = PetscOptionsBool("-pc_gamg_sym_graph",
184                             "Set for asymmetric matrices",
185                             "PCGAMGSetSymGraph",
186                             pc_gamg_agg->sym_graph,
187                             &pc_gamg_agg->sym_graph,
188                             &flag);CHKERRQ(ierr);
189 
190     /* -pc_gamg_square_graph */
191     pc_gamg_agg->square_graph = PETSC_TRUE;
192     ierr = PetscOptionsBool("-pc_gamg_square_graph",
193                             "For faster coarsening and lower coarse grid complexity",
194                             "PCGAMGSetSquareGraph",
195                             pc_gamg_agg->square_graph,
196                             &pc_gamg_agg->square_graph,
197                             &flag);CHKERRQ(ierr);
198   }
199   ierr = PetscOptionsTail();CHKERRQ(ierr);
200 
201   PetscFunctionReturn(0);
202 }
203 
204 /* -------------------------------------------------------------------------- */
205 /*
206    PCDestroy_AGG
207 
208   Input Parameter:
209    . pc -
210 */
211 #undef __FUNCT__
212 #define __FUNCT__ "PCDestroy_AGG"
213 PetscErrorCode PCDestroy_AGG(PC pc)
214 {
215   PetscErrorCode  ierr;
216   PC_MG           *mg = (PC_MG*)pc->data;
217   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
218   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
219 
220   PetscFunctionBegin;
221   if (pc_gamg_agg) {
222     ierr = PetscFree(pc_gamg_agg);CHKERRQ(ierr);
223     pc_gamg_agg = 0;
224   }
225 
226   /* call base class */
227   ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
228 
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;     assert(pc_gamg->data_cell_cols>0);
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       /* rotational modes */
291       if (coords) {
292         if (ndm == 2){
293           data += 2*M;
294           data[0] = -coords[2*kk+1];
295           data[1] =  coords[2*kk];
296         } else {
297           data += 3*M;
298           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
299           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
300           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
301         }
302       }
303     }
304   }
305 
306   pc_gamg->data_sz = arrsz;
307 
308   PetscFunctionReturn(0);
309 }
310 EXTERN_C_END
311 
312 typedef PetscInt NState;
313 static const NState NOT_DONE=-2;
314 static const NState DELETED=-1;
315 static const NState REMOVED=-3;
316 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
317 
318 /* -------------------------------------------------------------------------- */
319 /*
320    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
321      - AGG-MG specific: clears singletons out of 'selected_2'
322 
323    Input Parameter:
324    . Gmat_2 - glabal matrix of graph (data not defined)
325    . Gmat_1 - base graph to grab with
326    Input/Output Parameter:
327    . aggs_2 - linked list of aggs with gids)
328 */
329 #undef __FUNCT__
330 #define __FUNCT__ "smoothAggs"
331 static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */
332                                   const Mat Gmat_1, /* base graph */
333                                   /* const IS selected_2, [nselected local] selected vertices */
334                                   PetscCoarsenData *aggs_2 /* [nselected local] global ID of aggregate */
335                                  )
336 {
337   PetscErrorCode ierr;
338   PetscBool      isMPI;
339   Mat_SeqAIJ    *matA_1, *matB_1=0, *matA_2, *matB_2=0;
340   MPI_Comm       wcomm = ((PetscObject)Gmat_2)->comm;
341   PetscMPIInt    mype,npe;
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 = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
352   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
353   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
354 
355   if (PETSC_FALSE) {
356     PetscViewer viewer; char fname[32]; static int llev=0;
357     sprintf(fname,"Gmat2_%d.m",llev++);
358     PetscViewerASCIIOpen(wcomm,fname,&viewer);
359     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
360     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
361     ierr = PetscViewerDestroy(&viewer);
362   }
363 
364   /* get submatrices */
365   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
366   if (isMPI) {
367     /* grab matrix objects */
368     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
369     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
370     matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data;
371     matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data;
372     matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data;
373     matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data;
374 
375     /* force compressed row storage for B matrix in AuxMat */
376     matB_1->compressedrow.check = PETSC_TRUE;
377     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
378 
379     ierr = PetscMalloc(nloc*sizeof(PetscInt), &lid_cprowID_1);CHKERRQ(ierr);
380     for (lid = 0 ; lid < nloc ; lid++) lid_cprowID_1[lid] = -1;
381     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
382       PetscInt lid = matB_1->compressedrow.rindex[ix];
383       lid_cprowID_1[lid] = ix;
384     }
385   } else {
386     matA_1 = (Mat_SeqAIJ*)Gmat_1->data;
387     matA_2 = (Mat_SeqAIJ*)Gmat_2->data;
388     lid_cprowID_1 = PETSC_NULL;
389   }
390   assert(matA_1 && !matA_1->compressedrow.use);
391   assert(matB_1==0 || matB_1->compressedrow.use);
392   assert(matA_2 && !matA_2->compressedrow.use);
393   assert(matB_2==0 || matB_2->compressedrow.use);
394 
395   /* get state of locals and selected gid for deleted */
396   ierr = PetscMalloc(nloc*sizeof(NState), &lid_state);CHKERRQ(ierr);
397   ierr = PetscMalloc(nloc*sizeof(PetscScalar), &lid_parent_gid);CHKERRQ(ierr);
398   for (lid = 0 ; lid < nloc ; lid++) {
399     lid_parent_gid[lid] = -1.0;
400     lid_state[lid] = DELETED;
401   }
402 
403   /* set lid_state */
404   for (lid = 0 ; lid < nloc ; lid++) {
405     PetscCDPos pos;
406     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
407     if (pos) {
408       PetscInt gid1;
409       ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); assert(gid1==lid+my0);
410       lid_state[lid] = gid1;
411     }
412   }
413 
414   /* map local to selected local, DELETED means a ghost owns it */
415   for (lid=kk=0;lid<nloc;lid++){
416     NState state = lid_state[lid];
417     if (IS_SELECTED(state)){
418       PetscCDPos pos;
419       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
420       while(pos){
421         PetscInt gid1;
422         ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
423         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
424 
425         if (gid1 >= my0 && gid1 < Iend){
426           lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
427         }
428       }
429     }
430   }
431   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
432   if (isMPI) {
433     Vec          tempVec;
434     /* get 'cpcol_1_state' */
435     ierr = MatGetVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
436     for (kk=0,j=my0;kk<nloc;kk++,j++){
437       PetscScalar v = (PetscScalar)lid_state[kk];
438       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
439     }
440     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
441     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
442     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
443     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
444     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
445     /* get 'cpcol_2_state' */
446     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
447     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
448     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
449     /* get 'cpcol_2_par_orig' */
450     for (kk=0,j=my0;kk<nloc;kk++,j++){
451       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
452       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
453     }
454     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
455     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
456     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
457     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
458     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
459     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
460 
461     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
462   } /* ismpi */
463 
464   /* doit */
465   for (lid=0;lid<nloc;lid++){
466     NState state = lid_state[lid];
467     if (IS_SELECTED(state)) {
468       /* steal locals */
469       ii = matA_1->i; n = ii[lid+1] - ii[lid];
470       idx = matA_1->j + ii[lid];
471       for (j=0; j<n; j++) {
472         PetscInt lidj = idx[j], sgid;
473         NState statej = lid_state[lidj];
474         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
475           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
476           if (sgid >= my0 && sgid < Iend){       /* I'm stealing this local from a local sgid */
477             PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0;
478             PetscCDPos pos,last=PETSC_NULL;
479             /* looking for local from local so id_llist_2 works */
480             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
481             while(pos){
482               PetscInt gid;
483               ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr);
484               if (gid == gidj) {
485                 assert(last);
486                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
487                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
488                 hav = 1;
489                 break;
490               } else last = pos;
491 
492               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
493             }
494             if (hav!=1){
495               if (hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
496               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
497             }
498           } else{            /* I'm stealing this local, owned by a ghost */
499             if (sgid != -1) {
500               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.");
501             }
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=PETSC_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                   assert(last);
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=PETSC_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 
659   PetscFunctionReturn(0);
660 }
661 
662 /* -------------------------------------------------------------------------- */
663 /*
664    PCSetData_AGG - called if data is not set with PCSetCoordinates.
665       Looks in Mat for near null space.
666       Does not work for Stokes
667 
668   Input Parameter:
669    . pc -
670    . a_A - matrix to get (near) null space out of.
671 */
672 #undef __FUNCT__
673 #define __FUNCT__ "PCSetData_AGG"
674 PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
675 {
676   PetscErrorCode ierr;
677   PC_MG          *mg = (PC_MG*)pc->data;
678   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
679   MatNullSpace   mnull;
680 
681   PetscFunctionBegin;
682   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
683   if (!mnull) {
684     PetscInt bs,NN,MM;
685     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
686     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);  assert(MM%bs==0);
687     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, PETSC_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,PETSC_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     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); /* this does not work for Stokes */
705     pc_gamg->data_cell_rows = bs;
706   }
707   PetscFunctionReturn(0);
708 }
709 
710 /* -------------------------------------------------------------------------- */
711 /*
712  formProl0
713 
714    Input Parameter:
715    . agg_llists - list of arrays with aggregates
716    . bs - block size
717    . nSAvec - column bs of new P
718    . my0crs - global index of start of locals
719    . data_stride - bs*(nloc nodes + ghost nodes)
720    . data_in[data_stride*nSAvec] - local data on fine grid
721    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
722   Output Parameter:
723    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
724    . a_Prol - prolongation operator
725 */
726 #undef __FUNCT__
727 #define __FUNCT__ "formProl0"
728 static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists,/* list from selected vertices of aggregate unselected vertices */
729                                 const PetscInt bs,          /* (row) block size */
730                                 const PetscInt nSAvec,      /* column bs */
731                                 const PetscInt my0crs,      /* global index of start of locals */
732                                 const PetscInt data_stride, /* (nloc+nghost)*bs */
733                                 PetscReal      data_in[],   /* [data_stride][nSAvec] */
734                                 const PetscInt flid_fgid[], /* [data_stride/bs] */
735                                 PetscReal **a_data_out,
736                                 Mat a_Prol /* prolongation operator (output)*/
737                                )
738 {
739   PetscErrorCode ierr;
740   PetscInt  Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
741   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
742   PetscMPIInt    mype, npe;
743   PetscReal      *out_data;
744   PetscCDPos         pos;
745   GAMGHashTable  fgid_flid;
746 
747 /* #define OUT_AGGS */
748 #ifdef OUT_AGGS
749   static PetscInt llev = 0; char fname[32]; FILE *file = PETSC_NULL; PetscInt pM;
750 #endif
751 
752   PetscFunctionBegin;
753   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
754   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
755   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
756   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
757   Iend /= bs;
758   nghosts = data_stride/bs - nloc;
759 
760   ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
761   for (kk=0;kk<nghosts;kk++) {
762     ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
763   }
764 
765 #ifdef OUT_AGGS
766   sprintf(fname,"aggs_%d_%d.m",llev++,mype);
767   if (llev==1) {
768     file = fopen(fname,"w");
769   }
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   assert((ii/nSAvec)==my0crs); assert(nSelected==(jj-ii)/nSAvec);
781 
782   /* aloc space for coarse point data (output) */
783   out_data_stride = nSelected*nSAvec;
784   ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr);
785   for (ii=0;ii<out_data_stride*nSAvec;ii++) {
786     out_data[ii]=1.e300;
787   }
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           assert(flid>=0);
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 #ifdef 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*mype)%6]; str[13] = sim[(clid+17*mype)%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       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             assert(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       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",mype,__FUNCT__,ii,kk/bs,ndone,jj); */
906 
907 #ifdef OUT_AGGS
908   if (llev==1) fclose(file);
909 #endif
910   ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr);
911 
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,
928                                 const Mat Amat,
929                                 Mat *a_Gmat
930                                )
931 {
932   PetscErrorCode ierr;
933   PC_MG          *mg = (PC_MG*)pc->data;
934   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
935   const PetscInt verbose = pc_gamg->verbose;
936   const PetscReal vfilter = pc_gamg->threshold;
937   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
938   PetscMPIInt    mype,npe;
939   Mat            Gmat;
940   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
941   PetscBool  /* set,flg , */symm;
942 
943   PetscFunctionBegin;
944 #if defined PETSC_USE_LOG
945   ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
946 #endif
947   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
948   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
949 
950   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
951   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
952 
953   ierr  = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
954   ierr  = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
955 
956   *a_Gmat = Gmat;
957 
958 #if defined PETSC_USE_LOG
959   ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr);
960 #endif
961   PetscFunctionReturn(0);
962 }
963 
964 /* -------------------------------------------------------------------------- */
965 /*
966    PCGAMGCoarsen_AGG
967 
968   Input Parameter:
969    . a_pc - this
970   Input/Output Parameter:
971    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
972   Output Parameter:
973    . agg_lists - list of aggregates
974 */
975 #undef __FUNCT__
976 #define __FUNCT__ "PCGAMGCoarsen_AGG"
977 PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,
978                                   Mat *a_Gmat1,
979                                   PetscCoarsenData **agg_lists
980                                  )
981 {
982   PetscErrorCode ierr;
983   PC_MG          *mg = (PC_MG*)a_pc->data;
984   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
985   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
986   Mat             mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */
987   const PetscInt verbose = pc_gamg->verbose;
988   IS              perm;
989   PetscInt        Ii,nloc,bs,n,m;
990   PetscInt *permute;
991   PetscBool *bIndexSet;
992   MatCoarsen crs;
993   MPI_Comm        wcomm = ((PetscObject)Gmat1)->comm;
994   PetscMPIInt     mype,npe;
995 
996   PetscFunctionBegin;
997 #if defined PETSC_USE_LOG
998   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
999 #endif
1000   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
1001   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
1002   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
1003   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); assert(bs==1);
1004   nloc = n/bs;
1005 
1006   if (pc_gamg_agg->square_graph) {
1007     if (verbose > 1) {
1008       PetscPrintf(wcomm,"[%d]%s square graph\n",mype,__FUNCT__);
1009     }
1010     /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */
1011     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
1012     if (verbose > 2) {
1013       ierr = PetscPrintf(wcomm,"[%d]%s square graph done\n",mype,__FUNCT__);CHKERRQ(ierr);
1014       /* check for symetry */
1015       if (verbose > 4) {
1016 
1017       }
1018     }
1019   } else Gmat2 = Gmat1;
1020 
1021   /* get MIS aggs */
1022   /* randomize */
1023   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
1024   ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
1025   for (Ii = 0; Ii < nloc ; Ii++){
1026     bIndexSet[Ii] = PETSC_FALSE;
1027     permute[Ii] = Ii;
1028   }
1029   srand(1); /* make deterministic */
1030   for (Ii = 0; Ii < nloc ; Ii++) {
1031     PetscInt iSwapIndex = rand()%nloc;
1032     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1033       PetscInt iTemp = permute[iSwapIndex];
1034       permute[iSwapIndex] = permute[Ii];
1035       permute[Ii] = iTemp;
1036       bIndexSet[iSwapIndex] = PETSC_TRUE;
1037     }
1038   }
1039   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
1040 
1041   if (verbose > 1) {
1042     PetscPrintf(wcomm,"[%d]%s coarsen graph\n",mype,__FUNCT__);
1043   }
1044 
1045   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
1046 #if defined PETSC_GAMG_USE_LOG
1047   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1048 #endif
1049   ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
1050   /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */
1051   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1052   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1053   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1054   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
1055   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1056   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1057   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1058   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1059 
1060   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1061   ierr = PetscFree(permute);CHKERRQ(ierr);
1062 #if defined PETSC_GAMG_USE_LOG
1063   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1064 #endif
1065   if (verbose > 2) {
1066     PetscPrintf(wcomm,"[%d]%s coarsen graph done\n",mype,__FUNCT__);
1067   }
1068 
1069   /* smooth aggs */
1070   if (Gmat2 != Gmat1) {
1071     const PetscCoarsenData *llist = *agg_lists;
1072     ierr = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1073     ierr = MatDestroy(&Gmat1);CHKERRQ(ierr);
1074     *a_Gmat1 = Gmat2; /* output */
1075     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
1076     if (mat) SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1077   } else {
1078     const PetscCoarsenData *llist = *agg_lists;
1079     /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */
1080     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
1081     if (mat) {
1082       ierr = MatDestroy(&Gmat1);CHKERRQ(ierr);
1083       *a_Gmat1 = mat; /* output */
1084     }
1085   }
1086 #if defined PETSC_USE_LOG
1087   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1088 #endif
1089   if (verbose > 2) {
1090     PetscPrintf(wcomm,"[%d]%s PCGAMGCoarsen_AGG done\n",mype,__FUNCT__);
1091   }
1092 
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 /* -------------------------------------------------------------------------- */
1097 /*
1098  PCGAMGProlongator_AGG
1099 
1100  Input Parameter:
1101  . pc - this
1102  . Amat - matrix on this fine level
1103  . Graph - used to get ghost data for nodes in
1104  . agg_lists - list of aggregates
1105  Output Parameter:
1106  . a_P_out - prolongation operator to the next level
1107  */
1108 #undef __FUNCT__
1109 #define __FUNCT__ "PCGAMGProlongator_AGG"
1110 PetscErrorCode PCGAMGProlongator_AGG(PC pc,
1111                                       const Mat Amat,
1112                                       const Mat Gmat,
1113                                       PetscCoarsenData *agg_lists,
1114                                       Mat *a_P_out
1115                                      )
1116 {
1117   PC_MG          *mg = (PC_MG*)pc->data;
1118   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1119   const PetscInt verbose = pc_gamg->verbose;
1120   const PetscInt data_cols = pc_gamg->data_cell_cols;
1121   PetscErrorCode ierr;
1122   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1123   Mat            Prol;
1124   PetscMPIInt    mype, npe;
1125   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1126   const PetscInt col_bs = data_cols;
1127   PetscReal      *data_w_ghost;
1128   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1129 
1130   PetscFunctionBegin;
1131 #if defined PETSC_USE_LOG
1132   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1133 #endif
1134   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
1135   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
1136   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1137   ierr  = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
1138   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
1139 
1140   /* get 'nLocalSelected' */
1141   for (ii=0, nLocalSelected = 0 ; ii < nloc ; ii++){
1142     PetscBool ise;
1143     /* filter out singletons 0 or 1? */
1144     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1145     if (!ise) nLocalSelected++;
1146   }
1147 
1148   /* create prolongator, create P matrix */
1149   ierr = MatCreate(wcomm, &Prol);CHKERRQ(ierr);
1150   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1151   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1152   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
1153   ierr = MatSeqAIJSetPreallocation(Prol, data_cols, PETSC_NULL);CHKERRQ(ierr);
1154   ierr = MatMPIAIJSetPreallocation(Prol,data_cols, PETSC_NULL,data_cols, PETSC_NULL);CHKERRQ(ierr);
1155   /* nloc*bs, nLocalSelected*col_bs, */
1156   /* PETSC_DETERMINE, PETSC_DETERMINE, */
1157   /* data_cols, PETSC_NULL, data_cols, PETSC_NULL, */
1158   /* &Prol); */
1159 
1160   /* can get all points "removed" */
1161   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
1162   if (ii==0) {
1163     if (verbose > 0) {
1164       PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__);
1165     }
1166     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
1167     *a_P_out = PETSC_NULL;  /* out */
1168     PetscFunctionReturn(0);
1169   }
1170   if (verbose > 0) {
1171     PetscPrintf(wcomm,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/col_bs);
1172   }
1173   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
1174 
1175   assert((kk-myCrs0)%col_bs==0);
1176   myCrs0 = myCrs0/col_bs;
1177   assert((kk/col_bs-myCrs0)==nLocalSelected);
1178 
1179   /* create global vector of data in 'data_w_ghost' */
1180 #if defined PETSC_GAMG_USE_LOG
1181   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1182 #endif
1183   if (npe > 1) { /*  */
1184     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1185     ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr);
1186     for (jj = 0 ; jj < data_cols ; jj++){
1187       for (kk = 0 ; kk < bs ; kk++) {
1188         PetscInt ii,stride;
1189         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
1190         for (ii = 0 ; ii < nloc ; ii++, tp += bs){
1191           tmp_ldata[ii] = *tp;
1192         }
1193         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1194 
1195         if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1196           ierr = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr);
1197           nbnodes = bs*stride;
1198         }
1199         tp2 = data_w_ghost + jj*bs*stride + kk;
1200         for (ii = 0 ; ii < stride ; ii++, tp2 += bs){
1201           *tp2 = tmp_gdata[ii];
1202         }
1203         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
1204       }
1205     }
1206     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1207   } else {
1208     nbnodes = bs*nloc;
1209     data_w_ghost = (PetscReal*)pc_gamg->data;
1210   }
1211 
1212   /* get P0 */
1213   if (npe > 1){
1214     PetscReal *fid_glid_loc,*fiddata;
1215     PetscInt stride;
1216 
1217     ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr);
1218     for (kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1219     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1220     ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1221     for (kk=0;kk<stride;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
1222     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1223 
1224     assert(stride==nbnodes/bs);
1225     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1226   } else {
1227     ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr);
1228     for (kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk;
1229   }
1230 #if defined PETSC_GAMG_USE_LOG
1231   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1232   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1233 #endif
1234   {
1235     PetscReal *data_out = PETSC_NULL;
1236     ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1237     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
1238     pc_gamg->data = data_out;
1239     pc_gamg->data_cell_rows = data_cols;
1240     pc_gamg->data_sz = data_cols*data_cols*nLocalSelected;
1241   }
1242 #if defined PETSC_GAMG_USE_LOG
1243   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1244 #endif
1245   if (npe > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);
1246   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
1247 
1248   *a_P_out = Prol;  /* out */
1249 #if defined PETSC_USE_LOG
1250   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1251 #endif
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 /* -------------------------------------------------------------------------- */
1256 /*
1257    PCGAMGOptprol_AGG
1258 
1259   Input Parameter:
1260    . pc - this
1261    . Amat - matrix on this fine level
1262  In/Output Parameter:
1263    . a_P_out - prolongation operator to the next level
1264 */
1265 #undef __FUNCT__
1266 #define __FUNCT__ "PCGAMGOptprol_AGG"
1267 PetscErrorCode PCGAMGOptprol_AGG(PC pc,
1268                                   const Mat Amat,
1269                                   Mat *a_P
1270                                  )
1271 {
1272   PetscErrorCode ierr;
1273   PC_MG          *mg = (PC_MG*)pc->data;
1274   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1275   const PetscInt verbose = pc_gamg->verbose;
1276   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1277   PetscInt       jj;
1278   PetscMPIInt    mype,npe;
1279   Mat            Prol = *a_P;
1280   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1281 
1282   PetscFunctionBegin;
1283 #if defined PETSC_USE_LOG
1284   ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
1285 #endif
1286   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
1287   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
1288 
1289   /* smooth P0 */
1290   for (jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++){
1291     Mat tMat;
1292     Vec diag;
1293     PetscReal alpha, emax, emin;
1294 #if defined PETSC_GAMG_USE_LOG
1295     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1296 #endif
1297     if (jj == 0) {
1298       KSP eksp;
1299       Vec bb, xx;
1300       PC pc;
1301       ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr);
1302       ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr);
1303       {
1304         PetscRandom    rctx;
1305         ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
1306         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1307         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
1308         ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1309       }
1310 
1311       /* zeroing out BC rows -- needed for crazy matrices */
1312       {
1313         PetscInt Istart,Iend,ncols,jj,Ii;
1314         PetscScalar zero = 0.0;
1315         ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1316         for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) {
1317           ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1318           if(ncols <= 1) {
1319             ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
1320           }
1321           ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr);
1322         }
1323         ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
1324         ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
1325       }
1326 
1327       ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr);
1328       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
1329       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1330       ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1331       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1332 
1333       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
1334       ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1335       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
1336 
1337       ierr = KSPGetPC(eksp, &pc);CHKERRQ(ierr);
1338       ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1339 
1340       /* solve - keep stuff out of logging */
1341       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
1342       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
1343       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
1344       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
1345       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
1346 
1347       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1348       if (verbose > 0) {
1349         PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
1350                     __FUNCT__,emax,emin,PCJACOBI);
1351       }
1352       ierr = VecDestroy(&xx);CHKERRQ(ierr);
1353       ierr = VecDestroy(&bb);CHKERRQ(ierr);
1354       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
1355 
1356       if (pc_gamg->emax_id == -1) {
1357         ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr);
1358         assert(pc_gamg->emax_id != -1);
1359       }
1360       ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr);
1361     }
1362 
1363     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
1364     ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
1365     ierr = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr);
1366     ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
1367     ierr = VecReciprocal(diag);CHKERRQ(ierr);
1368     ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
1369     ierr = VecDestroy(&diag);CHKERRQ(ierr);
1370     alpha = -1.4/emax;
1371     ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
1372     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
1373     Prol = tMat;
1374 #if defined PETSC_GAMG_USE_LOG
1375     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1376 #endif
1377   }
1378 #if defined PETSC_USE_LOG
1379   ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr);
1380 #endif
1381   *a_P = Prol;
1382 
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 /* -------------------------------------------------------------------------- */
1387 /*
1388    PCGAMGKKTProl_AGG
1389 
1390   Input Parameter:
1391    . pc - this
1392    . Prol11 - matrix on this fine level
1393    . A21 - matrix on this fine level
1394  In/Output Parameter:
1395    . a_P22 - prolongation operator to the next level
1396 */
1397 #undef __FUNCT__
1398 #define __FUNCT__ "PCGAMGKKTProl_AGG"
1399 PetscErrorCode PCGAMGKKTProl_AGG(PC pc,
1400                                   const Mat Prol11,
1401                                   const Mat A21,
1402                                   Mat *a_P22
1403                                  )
1404 {
1405   PetscErrorCode ierr;
1406   PC_MG          *mg = (PC_MG*)pc->data;
1407   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1408   const PetscInt verbose = pc_gamg->verbose;
1409   /* PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;  */
1410   PetscMPIInt    mype,npe;
1411   Mat            Prol22,Tmat,Gmat;
1412   MPI_Comm       wcomm = ((PetscObject)pc)->comm;
1413   PetscCoarsenData *agg_lists;
1414 
1415   PetscFunctionBegin;
1416 #if defined PETSC_USE_LOG
1417   ierr = PetscLogEventBegin(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1418 #endif
1419   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
1420   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
1421 
1422   /* form C graph */
1423   ierr = MatMatMult(A21, Prol11, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Tmat);CHKERRQ(ierr);
1424   ierr = MatMatTransposeMult(Tmat, Tmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat);CHKERRQ(ierr);
1425   ierr = MatDestroy(&Tmat);CHKERRQ(ierr);
1426   ierr  = PCGAMGFilterGraph(&Gmat, 0.0, PETSC_FALSE, verbose);CHKERRQ(ierr);
1427 
1428   /* coarsen constraints */
1429   {
1430     MatCoarsen crs;
1431     ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
1432     ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
1433     ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
1434     ierr = MatCoarsenSetVerbose(crs, verbose);CHKERRQ(ierr);
1435     ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1436     ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1437     ierr = MatCoarsenGetData(crs, &agg_lists);CHKERRQ(ierr);
1438     ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1439   }
1440 
1441   /* form simple prolongation 'Prol22' */
1442   {
1443     PetscInt ii,mm,clid,my0,nloc,nLocalSelected;
1444     PetscScalar val = 1.0;
1445     /* get 'nLocalSelected' */
1446     ierr = MatGetLocalSize(Gmat, &nloc, &ii);CHKERRQ(ierr);
1447     for (ii=0, nLocalSelected = 0 ; ii < nloc ; ii++){
1448       PetscBool ise;
1449       /* filter out singletons 0 or 1? */
1450       ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1451       if (!ise) nLocalSelected++;
1452     }
1453 
1454     ierr = MatCreate(wcomm,&Prol22);CHKERRQ(ierr);
1455     ierr = MatSetSizes(Prol22,nloc, nLocalSelected,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1456     ierr = MatSetType(Prol22, MATAIJ);CHKERRQ(ierr);
1457     ierr = MatSeqAIJSetPreallocation(Prol22,1,PETSC_NULL);CHKERRQ(ierr);
1458     ierr = MatMPIAIJSetPreallocation(Prol22,1,PETSC_NULL,1,PETSC_NULL);CHKERRQ(ierr);
1459     /* ierr = MatCreateAIJ(wcomm, */
1460     /*                      nloc, nLocalSelected, */
1461     /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
1462     /*                      1, PETSC_NULL, 1, PETSC_NULL, */
1463     /*                      &Prol22); */
1464 
1465     ierr = MatGetOwnershipRange(Prol22, &my0, &ii);CHKERRQ(ierr);
1466     nloc = ii - my0;
1467 
1468     /* make aggregates */
1469     for (mm = clid = 0 ; mm < nloc ; mm++){
1470       ierr = PetscCDSizeAt(agg_lists, mm, &ii);CHKERRQ(ierr);
1471       if (ii > 0) {
1472         PetscInt asz=ii,cgid=my0+clid,rids[1000];
1473         PetscCDPos pos;
1474         if (asz>1000)SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Very large aggregate: %d",asz);
1475         ii = 0;
1476         ierr = PetscCDGetHeadPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1477         while(pos){
1478           PetscInt gid1;
1479           ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr);
1480           ierr = PetscCDGetNextPos(agg_lists,mm,&pos);CHKERRQ(ierr);
1481 
1482           rids[ii++] = gid1;
1483         }
1484         assert(ii==asz);
1485         /* add diagonal block of P0 */
1486         ierr = MatSetValues(Prol22,asz,rids,1,&cgid,&val,INSERT_VALUES);CHKERRQ(ierr);
1487 
1488         clid++;
1489       } /* coarse agg */
1490     } /* for all fine nodes */
1491     ierr = MatAssemblyBegin(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1492     ierr = MatAssemblyEnd(Prol22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1493   }
1494 
1495   /* clean up */
1496   ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
1497   ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
1498 #if defined PETSC_USE_LOG
1499   ierr = PetscLogEventEnd(PC_GAMGKKTProl_AGG,0,0,0,0);CHKERRQ(ierr);
1500 #endif
1501   *a_P22 = Prol22;
1502 
1503   PetscFunctionReturn(0);
1504 }
1505 
1506 /* -------------------------------------------------------------------------- */
1507 /*
1508    PCCreateGAMG_AGG
1509 
1510   Input Parameter:
1511    . pc -
1512 */
1513 #undef __FUNCT__
1514 #define __FUNCT__ "PCCreateGAMG_AGG"
1515 PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1516 {
1517   PetscErrorCode  ierr;
1518   PC_MG           *mg = (PC_MG*)pc->data;
1519   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1520   PC_GAMG_AGG      *pc_gamg_agg;
1521 
1522   PetscFunctionBegin;
1523 
1524   if(pc_gamg->subctx){
1525     /* call base class */
1526     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
1527   }
1528 
1529   /* create sub context for SA */
1530   ierr = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr);
1531   pc_gamg->subctx = pc_gamg_agg;
1532 
1533   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1534   pc->ops->destroy        = PCDestroy_AGG;
1535   /* reset does not do anything; setup not virtual */
1536 
1537   /* set internal function pointers */
1538   pc_gamg->graph = PCGAMGgraph_AGG;
1539   pc_gamg->coarsen = PCGAMGCoarsen_AGG;
1540   pc_gamg->prolongator = PCGAMGProlongator_AGG;
1541   pc_gamg->optprol = PCGAMGOptprol_AGG;
1542   pc_gamg->formkktprol = PCGAMGKKTProl_AGG;
1543 
1544   pc_gamg->createdefaultdata = PCSetData_AGG;
1545 
1546   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_AGG",PCSetCoordinates_AGG);CHKERRQ(ierr);
1547   PetscFunctionReturn(0);
1548 }
1549