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