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