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