xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision b43b03e97385fe8ad94dfdfbf5044a665bcc78db)
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 <private/kspimpl.h>
7 
8 #include <assert.h>
9 #include <petscblaslapack.h>
10 
11 typedef struct {
12   PetscInt nsmooths;
13   Mat aux_mat;
14   PetscBool sym_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 /* -------------------------------------------------------------------------- */
106 /*
107    PCSetFromOptions_GAMG_AGG
108 
109   Input Parameter:
110    . pc -
111 */
112 #undef __FUNCT__
113 #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
114 PetscErrorCode PCSetFromOptions_GAMG_AGG( PC pc )
115 {
116   PetscErrorCode  ierr;
117   PC_MG           *mg = (PC_MG*)pc->data;
118   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
119   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
120   PetscBool        flag;
121 
122   PetscFunctionBegin;
123   /* call base class */
124   ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr);
125 
126   ierr = PetscOptionsHead("GAMG-AGG options"); CHKERRQ(ierr);
127   {
128     /* -pc_gamg_agg_nsmooths */
129     pc_gamg_agg->nsmooths = 0;
130     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths",
131                            "smoothing steps for smoothed aggregation, usually 1 (0)",
132                            "PCGAMGSetNSmooths",
133                            pc_gamg_agg->nsmooths,
134                            &pc_gamg_agg->nsmooths,
135                            &flag);
136     CHKERRQ(ierr);
137 
138     /* -pc_gamg_sym_graph */
139     pc_gamg_agg->sym_graph = PETSC_FALSE;
140     ierr = PetscOptionsBool("-pc_gamg_sym_graph",
141                             "Set for asymetric matrices",
142                             "PCGAMGSetSymGraph",
143                             pc_gamg_agg->sym_graph,
144                             &pc_gamg_agg->sym_graph,
145                             &flag);
146     CHKERRQ(ierr);
147   }
148   ierr = PetscOptionsTail();CHKERRQ(ierr);
149 
150   if( pc_gamg->verbose > 1 ) {
151     PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done\n",0,__FUNCT__);
152   }
153 
154   PetscFunctionReturn(0);
155 }
156 
157 /* -------------------------------------------------------------------------- */
158 /*
159    PCDestroy_AGG
160 
161   Input Parameter:
162    . pc -
163 */
164 #undef __FUNCT__
165 #define __FUNCT__ "PCDestroy_AGG"
166 PetscErrorCode PCDestroy_AGG( PC pc )
167 {
168   PetscErrorCode  ierr;
169   PC_MG           *mg = (PC_MG*)pc->data;
170   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
171   PC_GAMG_AGG      *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
172 
173   PetscFunctionBegin;
174   if( pc_gamg_agg ) {
175     ierr = PetscFree(pc_gamg_agg);CHKERRQ(ierr);
176     pc_gamg_agg = 0;
177   }
178 
179   /* call base class */
180   ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr);
181 
182   PetscFunctionReturn(0);
183 }
184 
185 /* -------------------------------------------------------------------------- */
186 /*
187    PCSetCoordinates_AGG
188 
189    Input Parameter:
190    .  pc - the preconditioner context
191 */
192 EXTERN_C_BEGIN
193 #undef __FUNCT__
194 #define __FUNCT__ "PCSetCoordinates_AGG"
195 PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscReal *coords )
196 {
197   PC_MG          *mg = (PC_MG*)pc->data;
198   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
199   PetscErrorCode ierr;
200   PetscInt       arrsz,bs,my0,kk,ii,jj,nloc,Iend;
201   Mat            Amat = pc->pmat;
202 
203   PetscFunctionBegin;
204   PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 );
205   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
206   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
207   nloc = (Iend-my0)/bs;
208   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
209 
210   /* SA: null space vectors */
211   if( coords && bs==1 ) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
212   else if( coords ) pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* elasticity */
213   else pc_gamg->data_cell_cols = bs; /* no data, force SA with constant null space vectors */
214   pc_gamg->data_cell_rows = bs;
215 
216   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
217 
218   /* create data - syntactic sugar that should be refactored at some point */
219   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
220     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
221     ierr = PetscMalloc(arrsz*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
222   }
223   /* copy data in - column oriented */
224   for(kk=0;kk<nloc;kk++){
225     const PetscInt M = Iend - my0;
226     PetscReal *data = &pc_gamg->data[kk*bs];
227     if( pc_gamg->data_cell_cols==1 ) *data = 1.0;
228     else {
229       for(ii=0;ii<bs;ii++)
230         for(jj=0;jj<bs;jj++)
231           if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */
232           else data[ii*M + jj] = 0.0;
233       if( coords ) {
234         if( ndm == 2 ){ /* rotational modes */
235           data += 2*M;
236           data[0] = -coords[2*kk+1];
237           data[1] =  coords[2*kk];
238         }
239         else {
240           data += 3*M;
241           data[0] = 0.0;               data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
242           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;               data[2*M+1] =  coords[3*kk];
243           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
244         }
245       }
246     }
247   }
248 
249   pc_gamg->data_sz = arrsz;
250 
251   PetscFunctionReturn(0);
252 }
253 EXTERN_C_END
254 
255 typedef PetscInt NState;
256 static const NState NOT_DONE=-2;
257 static const NState DELETED=-1;
258 static const NState REMOVED=-3;
259 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
260 
261 /* -------------------------------------------------------------------------- */
262 /*
263    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
264      - AGG-MG specific: clears singletons out of 'selected_2'
265 
266    Input Parameter:
267    . Gmat_2 - glabal matrix of graph (data not defined)
268    . Gmat_1 - base graph to grab with
269    . selected_2 -
270    Input/Output Parameter:
271    . llist_aggs_2 - linked list of aggs, ghost lids are based on Gmat_2 (squared graph)
272 */
273 #undef __FUNCT__
274 #define __FUNCT__ "smoothAggs"
275 PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */
276                            const Mat Gmat_1, /* base graph, could be unsymmetic */
277                            const IS selected_2, /* [nselected total] selected vertices */
278                            IS llist_aggs_2 /* [nloc_nghost] global ID of aggregate */
279                            )
280 {
281   PetscErrorCode ierr;
282   PetscBool      isMPI;
283   Mat_SeqAIJ    *matA_1, *matB_1=0, *matA_2, *matB_2=0;
284   MPI_Comm       wcomm = ((PetscObject)Gmat_2)->comm;
285   PetscMPIInt    mype;
286   PetscInt       lid,*ii,*idx,ix,Iend,my0,nnodes_2,kk,n,j;
287   Mat_MPIAIJ    *mpimat_2 = 0, *mpimat_1=0;
288   const PetscInt nloc = Gmat_2->rmap->n;
289   PetscScalar   *cpcol_1_state,*cpcol_2_state,*deleted_parent_gid;
290   PetscInt      *lid_cprowID_1,*id_llist_2,*lid_cprowID_2;
291   NState        *lid_state;
292 
293   PetscFunctionBegin;
294   ierr = MPI_Comm_rank( wcomm, &mype );   CHKERRQ(ierr);
295   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);  CHKERRQ(ierr);
296 
297   if( !PETSC_TRUE ) {
298     PetscViewer viewer; char fname[32]; static int llev=0;
299     sprintf(fname,"Gmat2_%d.m",llev++);
300     PetscViewerASCIIOpen(wcomm,fname,&viewer);
301     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
302     ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr);
303     ierr = PetscViewerDestroy( &viewer );
304   }
305 
306   { /* copy linked list into temp buffer - should not work directly on pointer */
307     const PetscInt *llist_idx;
308     ierr = ISGetSize( llist_aggs_2, &nnodes_2 );        CHKERRQ(ierr);
309     ierr = PetscMalloc( nnodes_2*sizeof(PetscInt), &id_llist_2 ); CHKERRQ(ierr);
310     ierr = ISGetIndices( llist_aggs_2, &llist_idx );       CHKERRQ(ierr);
311     for(lid=0;lid<nnodes_2;lid++) id_llist_2[lid] = llist_idx[lid];
312     ierr = ISRestoreIndices( llist_aggs_2, &llist_idx );     CHKERRQ(ierr);
313   }
314 
315   /* get submatrices */
316   ierr = PetscTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr);
317   if(isMPI) {
318     /* grab matrix objects */
319     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
320     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
321     matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data;
322     matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data;
323     matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data;
324     matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data;
325 
326     /* force compressed row storage for B matrix in AuxMat */
327     matB_1->compressedrow.check = PETSC_TRUE;
328     ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);
329     CHKERRQ(ierr);
330 
331     ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_2 ); CHKERRQ(ierr);
332     ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr);
333     for(lid=0;lid<nloc;lid++) lid_cprowID_1[lid] = lid_cprowID_2[lid] = -1;
334     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
335       PetscInt lid = matB_1->compressedrow.rindex[ix];
336       lid_cprowID_1[lid] = ix;
337     }
338     for (ix=0; ix<matB_2->compressedrow.nrows; ix++) {
339       PetscInt lid = matB_2->compressedrow.rindex[ix];
340       lid_cprowID_2[lid] = ix;
341     }
342   }
343   else {
344     matA_1 = (Mat_SeqAIJ*)Gmat_1->data;
345     matA_2 = (Mat_SeqAIJ*)Gmat_2->data;
346     lid_cprowID_2 = lid_cprowID_1 = 0;
347   }
348   assert( matA_1 && !matA_1->compressedrow.use );
349   assert( matB_1==0 || matB_1->compressedrow.use );
350   assert( matA_2 && !matA_2->compressedrow.use );
351   assert( matB_2==0 || matB_2->compressedrow.use );
352 
353   /* get state of locals and selected gid for deleted */
354   ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr);
355   ierr = PetscMalloc( nloc*sizeof(PetscScalar), &deleted_parent_gid ); CHKERRQ(ierr);
356   for( lid = 0 ; lid < nloc ; lid++ ) {
357     deleted_parent_gid[lid] = -1.0;
358     lid_state[lid] = DELETED;
359   }
360   /* set index into compressed row 'lid_cprowID', not -1 means its a boundary node */
361   {
362     PetscInt nSelected;
363     const PetscInt *selected_idx;
364     /* set local selected */
365     ierr = ISGetSize( selected_2, &nSelected );        CHKERRQ(ierr);
366     ierr = ISGetIndices( selected_2, &selected_idx );       CHKERRQ(ierr);
367     for(kk=0;kk<nSelected;kk++){
368       PetscInt lid = selected_idx[kk];
369       if(lid<nloc) lid_state[lid] = (NState)(lid+my0); /* selected flag */
370       else break;
371       /* remove singletons */
372       ii = matA_2->i; n = ii[lid+1] - ii[lid];
373       if( n < 2 ) {
374         if(!lid_cprowID_2 || (ix=lid_cprowID_2[lid])==-1 || (matB_2->compressedrow.i[ix+1]-matB_2->compressedrow.i[ix])==0){
375           lid_state[lid] = REMOVED;
376         }
377       }
378     }
379     ierr = ISRestoreIndices( selected_2, &selected_idx );     CHKERRQ(ierr);
380   }
381   /* map local to selected local, -1 means a ghost owns it */
382   for(lid=kk=0;lid<nloc;lid++){
383     NState state = lid_state[lid];
384     if( IS_SELECTED(state) ){
385       PetscInt flid = lid;
386       do{
387         if(flid<nloc){
388           deleted_parent_gid[flid] = (PetscScalar)(lid + my0);
389         }
390         kk++;
391       } while( (flid=id_llist_2[flid]) != -1 );
392     }
393   }
394   /* get 'cpcol_1_state', 'cpcol_2_state' - uses mpimat_1->lvec & mpimat_2->lvec for temp space */
395   if (isMPI) {
396     Vec          tempVec;
397 
398     /* get 'cpcol_1_state' */
399     ierr = MatGetVecs( Gmat_1, &tempVec, 0 );         CHKERRQ(ierr);
400     for(kk=0,j=my0;kk<nloc;kk++,j++){
401       PetscScalar v = (PetscScalar)lid_state[kk];
402       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
403     }
404     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
405     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
406     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
407     CHKERRQ(ierr);
408     ierr =   VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);
409     CHKERRQ(ierr);
410     ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
411     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
412 
413     /* get 'cpcol_2_state' */
414     ierr = MatGetVecs( Gmat_2, &tempVec, 0 );         CHKERRQ(ierr);
415     for(kk=0,j=my0;kk<nloc;kk++,j++){
416       PetscScalar v = (PetscScalar)lid_state[kk];
417       ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES );  CHKERRQ(ierr);
418     }
419     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
420     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
421     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
422     CHKERRQ(ierr);
423     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
424     CHKERRQ(ierr);
425     ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
426     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
427   } /* ismpi */
428 
429   /* doit */
430   for(lid=0;lid<nloc;lid++){
431     NState state = lid_state[lid];
432     if( IS_SELECTED(state) ) {      /* steal locals */
433       ii = matA_1->i; n = ii[lid+1] - ii[lid];
434       idx = matA_1->j + ii[lid];
435       for (j=0; j<n; j++) {
436         PetscInt flid, lidj = idx[j], sgid;
437         NState statej = lid_state[lidj];
438         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(deleted_parent_gid[lidj])) != lid+my0) { /* steal local */
439           deleted_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this with _2 */
440           if( sgid >= my0 && sgid < my0+nloc ){       /* I'm stealing this local from a local */
441             PetscInt hav=0, flid2=sgid-my0, lastid;
442             /* looking for local from local so id_llist_2 works */
443             for( lastid=flid2, flid=id_llist_2[flid2] ; flid!=-1 ; flid=id_llist_2[flid] ) {
444               if( flid == lidj ) {
445                 id_llist_2[lastid] = id_llist_2[flid];                    /* remove lidj from list */
446                 id_llist_2[flid] = id_llist_2[lid]; id_llist_2[lid] = flid; /* insert 'lidj' into head of llist */
447                 hav++;
448                 break;
449               }
450               lastid = flid;
451             }
452             if(hav!=1){
453               if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
454               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav);
455             }
456           }
457           else{            /* I'm stealing this local, owned by a ghost - ok to use _2, local */
458             assert(sgid==-1);
459             id_llist_2[lidj] = id_llist_2[lid]; id_llist_2[lid] = lidj; /* insert 'lidj' into head of llist */
460             /* local remove at end, off add/rm at end */
461           }
462         }
463       }
464     }
465     else if( state == DELETED && lid_cprowID_1 ) {
466       PetscInt sgidold = (PetscInt)PetscRealPart(deleted_parent_gid[lid]);
467       /* see if I have a selected ghost neighbor that will steal me */
468       if( (ix=lid_cprowID_1[lid]) != -1 ){
469         ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
470         idx = matB_1->j + ii[ix];
471         for( j=0 ; j<n ; j++ ) {
472           PetscInt cpid = idx[j];
473           NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
474           if( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */
475             deleted_parent_gid[lid] = (PetscScalar)statej; /* send who selected with _2 */
476             if( sgidold>=my0 && sgidold<(my0+nloc) ) { /* this was mine */
477               PetscInt lastid,hav=0,flid,oldslidj=sgidold-my0;
478               /* remove from 'oldslidj' list, local so _2 is OK */
479               for( lastid=oldslidj, flid=id_llist_2[oldslidj] ; flid != -1 ; flid=id_llist_2[flid] ) {
480                 if( flid == lid ) {
481                   id_llist_2[lastid] = id_llist_2[flid];   /* remove lid from oldslidj list */
482                   hav++;
483                   break;
484                 }
485                 lastid = flid;
486               }
487               if(hav!=1){
488                 if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
489                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav);
490               }
491               id_llist_2[lid] = -1; /* terminate linked list - needed? */
492             }
493             else assert(id_llist_2[lid] == -1);
494           }
495         }
496       }
497     } /* selected/deleted */
498     else assert(state == REMOVED || !lid_cprowID_1);
499   } /* node loop */
500 
501   if( isMPI ) {
502     PetscScalar *cpcol_2_sel_gid;
503     Vec          tempVec;
504     PetscInt     cpid;
505 
506     ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr);
507     ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr);
508 
509     /* get 'cpcol_2_sel_gid' */
510     ierr = MatGetVecs( Gmat_2, &tempVec, 0 );         CHKERRQ(ierr);
511     for(kk=0,j=my0;kk<nloc;kk++,j++){
512       ierr = VecSetValues( tempVec, 1, &j, &deleted_parent_gid[kk], INSERT_VALUES );  CHKERRQ(ierr);
513     }
514     ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr);
515     ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr);
516     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
517     CHKERRQ(ierr);
518     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);
519     CHKERRQ(ierr);
520     ierr = VecDestroy( &tempVec ); CHKERRQ(ierr);
521 
522     ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr);
523 
524     /* look for deleted ghosts and see if they moved */
525     for(lid=0;lid<nloc;lid++){
526       NState state = lid_state[lid];
527       if( IS_SELECTED(state) ){
528         PetscInt flid,lastid,old_sgid=lid+my0;
529         /* look for deleted ghosts and see if they moved */
530         for( lastid=lid, flid=id_llist_2[lid] ; flid!=-1 ; flid=id_llist_2[flid] ) {
531           if( flid>=nloc ) {
532             PetscInt cpid = flid-nloc, sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]);
533             if( sgid_new != old_sgid && sgid_new != -1 ) {
534               id_llist_2[lastid] = id_llist_2[flid];                    /* remove 'flid' from list */
535               id_llist_2[flid] = -1;
536               flid = lastid;
537             } /* if it changed parents */
538             else lastid = flid;
539           } /* for ghost nodes */
540           else lastid = flid;
541         } /* loop over list of deleted */
542       } /* selected */
543     }
544 
545     /* look at ghosts, see if they changed, and moved here */
546     for(cpid=0;cpid<nnodes_2-nloc;cpid++){
547       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]);
548       if( sgid_new>=my0 && sgid_new<(my0+nloc) ) { /* this is mine */
549         PetscInt lastid,flid,slid_new=sgid_new-my0,flidj=nloc+cpid,hav=0;
550         for( lastid=slid_new, flid=id_llist_2[slid_new] ; flid != -1 ; flid=id_llist_2[flid] ) {
551           if( flid == flidj ) {
552             hav++;
553             break;
554           }
555           lastid = flid;
556         }
557         if( hav != 1 ){
558           assert(id_llist_2[flidj] == -1);
559           id_llist_2[flidj] = id_llist_2[slid_new]; id_llist_2[slid_new] = flidj; /* insert 'flidj' into head of llist */
560         }
561       }
562     }
563 
564     ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr);
565     ierr = PetscFree( lid_cprowID_1 );  CHKERRQ(ierr);
566     ierr = PetscFree( lid_cprowID_2 );  CHKERRQ(ierr);
567   }
568 
569   /* copy out new aggs */
570   ierr = ISGeneralSetIndices(llist_aggs_2, nnodes_2, id_llist_2, PETSC_COPY_VALUES ); CHKERRQ(ierr);
571 
572   ierr = PetscFree( id_llist_2 );  CHKERRQ(ierr);
573   ierr = PetscFree( deleted_parent_gid );  CHKERRQ(ierr);
574   ierr = PetscFree( lid_state );  CHKERRQ(ierr);
575 
576   PetscFunctionReturn(0);
577 }
578 
579 /* -------------------------------------------------------------------------- */
580 /*
581    PCSetData_AGG
582 
583   Input Parameter:
584    . pc -
585 */
586 #undef __FUNCT__
587 #define __FUNCT__ "PCSetData_AGG"
588 PetscErrorCode PCSetData_AGG( PC pc )
589 {
590   PetscErrorCode  ierr;
591   PetscFunctionBegin;
592   ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 
596 /* -------------------------------------------------------------------------- */
597 /*
598  formProl0
599 
600    Input Parameter:
601    . selected - list of selected local ID, includes selected ghosts
602    . locals_llist - linked list with aggregates
603    . bs - block size
604    . nSAvec - num columns of new P
605    . my0crs - global index of locals
606    . data_stride - bs*(nloc nodes + ghost nodes)
607    . data_in[data_stride*nSAvec] - local data on fine grid
608    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
609   Output Parameter:
610    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
611    . a_Prol - prolongation operator
612 */
613 #undef __FUNCT__
614 #define __FUNCT__ "formProl0"
615 PetscErrorCode formProl0( IS selected, /* list of selected local ID, includes selected ghosts */
616                           IS locals_llist, /* linked list from selected vertices of aggregate unselected vertices */
617                           const PetscInt bs,
618                           const PetscInt nSAvec,
619                           const PetscInt my0crs,
620                           const PetscInt data_stride,
621                           PetscReal data_in[],
622                           const PetscInt flid_fgid[],
623                           PetscReal **a_data_out,
624                           Mat a_Prol /* prolongation operator (output)*/
625                           )
626 {
627   PetscErrorCode ierr;
628   PetscInt  Istart,Iend,nFineLoc,clid,flid,aggID,kk,jj,ii,mm,nLocalSelected,ndone,nSelected,minsz;
629   MPI_Comm       wcomm = ((PetscObject)a_Prol)->comm;
630   PetscMPIInt    mype, npe;
631   const PetscInt *selected_idx,*llist_idx;
632   PetscReal      *out_data;
633 
634   PetscFunctionBegin;
635   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
636   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
637   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );    CHKERRQ(ierr);
638   nFineLoc = (Iend-Istart)/bs; assert((Iend-Istart)%bs==0);
639 
640   ierr = ISGetSize( selected, &nSelected );        CHKERRQ(ierr);
641   ierr = ISGetIndices( selected, &selected_idx );       CHKERRQ(ierr);
642   ierr = ISGetIndices( locals_llist, &llist_idx );      CHKERRQ(ierr);
643   for(kk=0,nLocalSelected=0;kk<nSelected;kk++){
644     PetscInt lid = selected_idx[kk];
645     if( lid<nFineLoc && llist_idx[lid] != -1 ) nLocalSelected++;
646   }
647 
648   /* aloc space for coarse point data (output) */
649 #define DATA_OUT_STRIDE (nLocalSelected*nSAvec)
650   ierr = PetscMalloc( DATA_OUT_STRIDE*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr);
651   for(ii=0;ii<DATA_OUT_STRIDE*nSAvec;ii++) out_data[ii]=1.e300;
652   *a_data_out = out_data; /* output - stride nLocalSelected*nSAvec */
653 
654   /* find points and set prolongation */
655   minsz = 100;
656   ndone = 0;
657   for( mm = clid = 0 ; mm < nSelected ; mm++ ){
658     PetscInt lid = selected_idx[mm];
659     PetscInt cgid = my0crs + clid, cids[100];
660 
661     if( lid>=nFineLoc || llist_idx[lid]==-1 ) continue; /* skip ghost or singleton */
662 
663     /* count agg */
664     aggID = 0;
665     flid = selected_idx[mm]; assert(flid != -1);
666     do{
667       aggID++;
668     } while( (flid=llist_idx[flid]) != -1 );
669     if( aggID<minsz ) minsz = aggID;
670 
671     /* get block */
672     {
673       PetscBLASInt   asz=aggID,M=asz*bs,N=nSAvec,INFO;
674       PetscBLASInt   Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs;
675       PetscScalar    *qqc,*qqr,*TAU,*WORK;
676       PetscInt       *fids;
677 
678       ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr);
679       ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr);
680       ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr);
681       ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr);
682       ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr);
683 
684       flid = selected_idx[mm];
685       aggID = 0;
686       do{
687         /* copy in B_i matrix - column oriented */
688         PetscReal *data = &data_in[flid*bs];
689         for( kk = ii = 0; ii < bs ; ii++ ) {
690           for( jj = 0; jj < N ; jj++ ) {
691             qqc[jj*Mdata + aggID*bs + ii] = data[jj*data_stride + ii];
692           }
693         }
694 
695         /* set fine IDs */
696         for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
697 
698         aggID++;
699       }while( (flid=llist_idx[flid]) != -1 );
700 
701       /* pad with zeros */
702       for( ii = asz*bs; ii < Mdata ; ii++ ) {
703 	for( jj = 0; jj < N ; jj++, kk++ ) {
704 	  qqc[jj*Mdata + ii] = .0;
705 	}
706       }
707 
708       ndone += aggID;
709       /* QR */
710       LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
711       if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error");
712       /* get R - column oriented - output B_{i+1} */
713       {
714         PetscReal *data = &out_data[clid*nSAvec];
715         for( jj = 0; jj < nSAvec ; jj++ ) {
716           for( ii = 0; ii < nSAvec ; ii++ ) {
717             assert(data[jj*DATA_OUT_STRIDE + ii] == 1.e300);
718             if( ii <= jj ) data[jj*DATA_OUT_STRIDE + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
719 	    else data[jj*DATA_OUT_STRIDE + ii] = 0.;
720           }
721         }
722       }
723 
724       /* get Q - row oriented */
725       LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO );
726       if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO);
727 
728       for( ii = 0 ; ii < M ; ii++ ){
729         for( jj = 0 ; jj < N ; jj++ ) {
730           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
731         }
732       }
733 
734       /* add diagonal block of P0 */
735       for(kk=0;kk<N;kk++) {
736         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
737       }
738       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr);
739 
740       ierr = PetscFree( qqc );  CHKERRQ(ierr);
741       ierr = PetscFree( qqr );  CHKERRQ(ierr);
742       ierr = PetscFree( TAU );  CHKERRQ(ierr);
743       ierr = PetscFree( WORK );  CHKERRQ(ierr);
744       ierr = PetscFree( fids );  CHKERRQ(ierr);
745     } /* scoping */
746     clid++;
747   } /* for all coarse nodes */
748 
749 /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */
750 /* MatGetSize( a_Prol, &kk, &jj ); */
751 /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */
752 /* PetscPrintf(PETSC_COMM_WORLD," **** [%d]%s %d total done, N=%d (%d local done), min agg. size = %d\n",mype,__FUNCT__,ii,kk/bs,ndone,jj); */
753 
754   ierr = ISRestoreIndices( selected, &selected_idx );     CHKERRQ(ierr);
755   ierr = ISRestoreIndices( locals_llist, &llist_idx );     CHKERRQ(ierr);
756   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
757   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758 
759   PetscFunctionReturn(0);
760 }
761 
762 /* -------------------------------------------------------------------------- */
763 /*
764    PCGAMGgraph_AGG
765 
766   Input Parameter:
767    . pc - this
768    . Amat - matrix on this fine level
769   Output Parameter:
770    . a_Gmat -
771 */
772 #undef __FUNCT__
773 #define __FUNCT__ "PCGAMGgraph_AGG"
774 PetscErrorCode PCGAMGgraph_AGG( PC pc,
775                                 const Mat Amat,
776                                 Mat *a_Gmat
777                                 )
778 {
779   PetscErrorCode ierr;
780   PC_MG          *mg = (PC_MG*)pc->data;
781   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
782   const PetscInt verbose = pc_gamg->verbose;
783   const PetscReal vfilter = pc_gamg->threshold;
784   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
785   PetscMPIInt    mype,npe;
786   Mat            Gmat, Gmat2;
787   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
788 
789   PetscFunctionBegin;
790   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
791   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
792 
793   ierr  = createSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr );
794 
795   ierr  = scaleFilterGraph( &Gmat, vfilter, pc_gamg_agg->sym_graph, verbose ); CHKERRQ( ierr );
796 
797   ierr = MatTransposeMatMult( Gmat, Gmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 );
798   CHKERRQ(ierr);
799 
800   /* attach auxilary matrix */
801   pc_gamg_agg->aux_mat = Gmat;
802 
803   *a_Gmat = Gmat2;
804 
805   PetscFunctionReturn(0);
806 }
807 
808 /* -------------------------------------------------------------------------- */
809 /*
810    PCGAMGCoarsen_AGG
811 
812   Input Parameter:
813    . pc - this
814    . Gmat2 - matrix on this fine level
815   Output Parameter:
816    . a_selected - prolongation operator to the next level
817    . a_llist_parent - data of coarse grid points (num local columns in 'a_P_out')
818 */
819 #undef __FUNCT__
820 #define __FUNCT__ "PCGAMGCoarsen_AGG"
821 PetscErrorCode PCGAMGCoarsen_AGG( PC pc,
822                                   const Mat Gmat2,
823                                   IS *a_selected,
824                                   IS *a_llist_parent
825                                   )
826 {
827   PetscErrorCode ierr;
828   PC_MG          *mg = (PC_MG*)pc->data;
829   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
830   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
831   Mat             Gmat1; /* unsquared graph (not symetrized!) */
832   IS              perm, selected, llist_parent;
833   PetscInt        Ii,nloc,bs,n,m;
834   PetscInt *permute;
835   PetscBool *bIndexSet;
836   MatCoarsen crs;
837   MPI_Comm        wcomm = ((PetscObject)Gmat2)->comm;
838   /* PetscMPIInt     mype,npe; */
839 
840   PetscFunctionBegin;
841   /* ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr); */
842   /* ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr); */
843   ierr = MatGetLocalSize( Gmat2, &n, &m ); CHKERRQ(ierr);
844   ierr = MatGetBlockSize( Gmat2, &bs ); CHKERRQ(ierr); assert(bs==1);
845   nloc = n/bs;
846 
847   /* get unsquared graph */
848   Gmat1 = pc_gamg_agg->aux_mat; pc_gamg_agg->aux_mat = 0;
849 
850   /* get MIS aggs */
851   /* randomize */
852   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
853   ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
854   for ( Ii = 0; Ii < nloc ; Ii++ ){
855     bIndexSet[Ii] = PETSC_FALSE;
856     permute[Ii] = Ii;
857   }
858   srand(1); /* make deterministic */
859   for ( Ii = 0; Ii < nloc ; Ii++ ) {
860     PetscInt iSwapIndex = rand()%nloc;
861     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
862       PetscInt iTemp = permute[iSwapIndex];
863       permute[iSwapIndex] = permute[Ii];
864       permute[Ii] = iTemp;
865       bIndexSet[iSwapIndex] = PETSC_TRUE;
866     }
867   }
868   ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
869 
870   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);
871   CHKERRQ(ierr);
872 #if defined PETSC_USE_LOG
873   ierr = PetscLogEventBegin(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
874 #endif
875 
876   //ierr = maxIndSetAgg( perm, Gmat2, PETSC_TRUE, pc_gamg->verbose, &selected, &llist_parent ); CHKERRQ(ierr);
877 
878   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
879   ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr);
880   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
881   ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr);
882   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
883   ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr);
884   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
885   ierr = MatCoarsenGetMISAggLists( crs, &selected, &llist_parent ); CHKERRQ(ierr);
886   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
887 
888   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
889   ierr = PetscFree( permute );  CHKERRQ(ierr);
890 #if defined PETSC_USE_LOG
891   ierr = PetscLogEventEnd(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
892 #endif
893   /* smooth aggs */
894   ierr = smoothAggs( Gmat2, Gmat1, selected, llist_parent ); CHKERRQ(ierr);
895 
896   ierr = MatDestroy( &Gmat1 );  CHKERRQ(ierr);
897 
898   *a_selected = selected;
899   *a_llist_parent = llist_parent;
900 
901   PetscFunctionReturn(0);
902 }
903 
904 /* -------------------------------------------------------------------------- */
905 /*
906  PCGAMGprolongator_AGG
907 
908  Input Parameter:
909  . pc - this
910  . Amat - matrix on this fine level
911  . Graph - used to get ghost data for nodes in
912  . selected - [nselected inc. chosts]
913  . llist_parent - [nloc + Gmat.nghost] linked list
914  Output Parameter:
915  . a_P_out - prolongation operator to the next level
916  */
917 #undef __FUNCT__
918 #define __FUNCT__ "PCGAMGprolongator_AGG"
919 PetscErrorCode PCGAMGprolongator_AGG( PC pc,
920                                      const Mat Amat,
921                                      const Mat Gmat,
922                                      IS selected,
923                                      IS llist_parent,
924                                      Mat *a_P_out
925                                      )
926 {
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 PetscInt data_cols = pc_gamg->data_cell_cols;
931   PetscErrorCode ierr;
932   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
933   Mat            Prol;
934   PetscMPIInt    mype, npe;
935   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
936   const PetscInt *selected_idx,*llist_idx,col_bs=data_cols;
937   PetscReal      *data_w_ghost;
938   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
939 
940   PetscFunctionBegin;
941   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
942   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
943   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
944   ierr  = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr );
945   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
946 
947   /* get 'nLocalSelected' */
948   ierr = ISGetSize( selected, &kk );        CHKERRQ(ierr);
949   ierr = ISGetIndices( selected, &selected_idx );     CHKERRQ(ierr);
950   ierr = ISGetIndices( llist_parent, &llist_idx );     CHKERRQ(ierr);
951   for(ii=0,nLocalSelected=0;ii<kk;ii++){
952     PetscInt lid = selected_idx[ii];
953     /* filter out singletons */
954     if( lid<nloc && llist_idx[lid] != -1) nLocalSelected++;
955   }
956   ierr = ISRestoreIndices( selected, &selected_idx );     CHKERRQ(ierr);
957   ierr = ISRestoreIndices( llist_parent, &llist_idx );     CHKERRQ(ierr);
958 
959   /* create prolongator, create P matrix */
960   ierr = MatCreateMPIAIJ( wcomm,
961                           nloc*bs, nLocalSelected*col_bs,
962                           PETSC_DETERMINE, PETSC_DETERMINE,
963                           data_cols, PETSC_NULL, data_cols, PETSC_NULL,
964                           &Prol );
965   CHKERRQ(ierr);
966 
967   /* can get all points "removed" */
968   ierr =  MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr);
969   if( ii==0 ) {
970     if( verbose ) {
971       PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__);
972     }
973     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
974     *a_P_out = PETSC_NULL;  /* out */
975     PetscFunctionReturn(0);
976   }
977   if( verbose ) {
978     PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/bs);
979   }
980   ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr);
981   myCrs0 = myCrs0/col_bs;
982 
983   /* create global vector of data in 'data_w_ghost' */
984 #if defined PETSC_USE_LOG
985   ierr = PetscLogEventBegin(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
986 #endif
987   if (npe > 1) { /*  */
988     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
989     ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr);
990     for( jj = 0 ; jj < data_cols ; jj++ ){
991       for( kk = 0 ; kk < bs ; kk++) {
992         PetscInt ii,nnodes;
993         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
994         for( ii = 0 ; ii < nloc ; ii++, tp += bs ){
995           tmp_ldata[ii] = *tp;
996         }
997         ierr = getDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata );
998         CHKERRQ(ierr);
999         if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */
1000           ierr = PetscMalloc( nnodes*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr);
1001           nbnodes = bs*nnodes;
1002         }
1003         tp2 = data_w_ghost + jj*bs*nnodes + kk;
1004         for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs ){
1005           *tp2 = tmp_gdata[ii];
1006         }
1007         ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr);
1008       }
1009     }
1010     ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr);
1011   }
1012   else {
1013     nbnodes = bs*nloc;
1014     data_w_ghost = (PetscReal*)pc_gamg->data;
1015   }
1016 
1017   /* get P0 */
1018   if( npe > 1 ){
1019     PetscReal *fid_glid_loc,*fiddata;
1020     PetscInt nnodes;
1021 
1022     ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr);
1023     for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1024     ierr = getDataWithGhosts( Gmat, 1, fid_glid_loc, &nnodes, &fiddata );
1025     CHKERRQ(ierr);
1026     ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
1027     for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
1028     ierr = PetscFree( fiddata ); CHKERRQ(ierr);
1029     assert(nnodes==nbnodes/bs);
1030     ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr);
1031   }
1032   else {
1033     ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr);
1034     for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk;
1035   }
1036 #if defined PETSC_USE_LOG
1037   ierr = PetscLogEventEnd(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
1038   ierr = PetscLogEventBegin(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1039 #endif
1040   {
1041     PetscReal *data_out;
1042     ierr = formProl0( selected, llist_parent, bs, data_cols, myCrs0, nbnodes,
1043                       data_w_ghost, flid_fgid, &data_out, Prol );
1044     CHKERRQ(ierr);
1045     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
1046     pc_gamg->data = data_out;
1047     pc_gamg->data_cell_rows = data_cols;
1048     pc_gamg->data_sz = data_cols*data_cols*nLocalSelected;
1049   }
1050  #if defined PETSC_USE_LOG
1051   ierr = PetscLogEventEnd(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1052 #endif
1053   if (npe > 1) ierr = PetscFree( data_w_ghost );      CHKERRQ(ierr);
1054   ierr = PetscFree( flid_fgid ); CHKERRQ(ierr);
1055 
1056   /* attach block size of columns */
1057   if( pc_gamg->col_bs_id == -1 ) {
1058     ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id ); assert(pc_gamg->col_bs_id != -1 );
1059   }
1060   ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr);
1061 
1062   *a_P_out = Prol;  /* out */
1063 
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 /* -------------------------------------------------------------------------- */
1068 /*
1069    PCGAMGoptprol_AGG
1070 
1071   Input Parameter:
1072    . pc - this
1073    . Amat - matrix on this fine level
1074  In/Output Parameter:
1075    . a_P_out - prolongation operator to the next level
1076 */
1077 #undef __FUNCT__
1078 #define __FUNCT__ "PCGAMGoptprol_AGG"
1079 PetscErrorCode PCGAMGoptprol_AGG( PC pc,
1080                                   const Mat Amat,
1081                                   Mat *a_P
1082                                   )
1083 {
1084   PetscErrorCode ierr;
1085   PC_MG          *mg = (PC_MG*)pc->data;
1086   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1087   const PetscInt verbose = pc_gamg->verbose;
1088   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1089   PetscInt       jj;
1090   PetscMPIInt    mype,npe;
1091   Mat            Prol = *a_P;
1092   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
1093 
1094   PetscFunctionBegin;
1095   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
1096   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
1097 
1098   /* smooth P0 */
1099   for( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){
1100     Mat tMat;
1101     Vec diag;
1102     PetscReal alpha, emax, emin;
1103 #if defined PETSC_USE_LOG
1104     ierr = PetscLogEventBegin(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1105 #endif
1106     if( jj == 0 ) {
1107       KSP eksp;
1108       Vec bb, xx;
1109       PC pc;
1110       ierr = MatGetVecs( Amat, &bb, 0 );         CHKERRQ(ierr);
1111       ierr = MatGetVecs( Amat, &xx, 0 );         CHKERRQ(ierr);
1112       {
1113         PetscRandom    rctx;
1114         ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
1115         ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1116         ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
1117         ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
1118       }
1119       ierr = KSPCreate(wcomm,&eksp);                            CHKERRQ(ierr);
1120       ierr = KSPAppendOptionsPrefix( eksp, "est_");         CHKERRQ(ierr);
1121       ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
1122       ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE );    CHKERRQ(ierr);
1123       ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN );
1124       CHKERRQ( ierr );
1125       ierr = KSPGetPC( eksp, &pc );                              CHKERRQ( ierr );
1126       ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr);  /* smoother */
1127       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);
1128       CHKERRQ(ierr);
1129       ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
1130       ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE );        CHKERRQ(ierr);
1131 
1132       /* solve - keep stuff out of logging */
1133       ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
1134       ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
1135       ierr = KSPSolve( eksp, bb, xx );                              CHKERRQ(ierr);
1136       ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
1137       ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
1138 
1139       ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
1140       if( verbose ) {
1141         PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n",
1142                     __FUNCT__,emax,emin,PCJACOBI);
1143       }
1144       ierr = VecDestroy( &xx );       CHKERRQ(ierr);
1145       ierr = VecDestroy( &bb );       CHKERRQ(ierr);
1146       ierr = KSPDestroy( &eksp );     CHKERRQ(ierr);
1147 
1148       if( pc_gamg->emax_id == -1 ) {
1149         ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id );
1150         assert(pc_gamg->emax_id != -1 );
1151       }
1152       ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr);
1153     }
1154 
1155     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
1156     ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat );   CHKERRQ(ierr);
1157     ierr = MatGetVecs( Amat, &diag, 0 );    CHKERRQ(ierr);
1158     ierr = MatGetDiagonal( Amat, diag );    CHKERRQ(ierr); /* effectively PCJACOBI */
1159     ierr = VecReciprocal( diag );         CHKERRQ(ierr);
1160     ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr);
1161     ierr = VecDestroy( &diag );           CHKERRQ(ierr);
1162     alpha = -1.5/emax;
1163     ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN );           CHKERRQ(ierr);
1164     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
1165     Prol = tMat;
1166 #if defined PETSC_USE_LOG
1167     ierr = PetscLogEventEnd(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
1168 #endif
1169   }
1170 
1171   *a_P = Prol;
1172 
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 /* -------------------------------------------------------------------------- */
1177 /*
1178    PCCreateGAMG_AGG
1179 
1180   Input Parameter:
1181    . pc -
1182 */
1183 #undef __FUNCT__
1184 #define __FUNCT__ "PCCreateGAMG_AGG"
1185 PetscErrorCode  PCCreateGAMG_AGG( PC pc )
1186 {
1187   PetscErrorCode  ierr;
1188   PC_MG           *mg = (PC_MG*)pc->data;
1189   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1190   PC_GAMG_AGG      *pc_gamg_agg;
1191 
1192   PetscFunctionBegin;
1193   /* create sub context for SA */
1194   ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr);
1195   assert(!pc_gamg->subctx);
1196   pc_gamg->subctx = pc_gamg_agg;
1197 
1198   pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1199   pc->ops->destroy        = PCDestroy_AGG;
1200   /* reset does not do anything; setup not virtual */
1201 
1202   /* set internal function pointers */
1203   pc_gamg->graph = PCGAMGgraph_AGG;
1204   pc_gamg->coarsen = PCGAMGCoarsen_AGG;
1205   pc_gamg->prolongator = PCGAMGprolongator_AGG;
1206   pc_gamg->optprol = PCGAMGoptprol_AGG;
1207 
1208   pc_gamg->createdefaultdata = PCSetData_AGG;
1209 
1210   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1211                                             "PCSetCoordinates_C",
1212                                             "PCSetCoordinates_AGG",
1213                                             PCSetCoordinates_AGG);
1214   PetscFunctionReturn(0);
1215 }
1216