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