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