xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 4e17fbf5bafa47185df4021e6c4f5fb738d34ecb)
1 /*
2  GAMG geometric-algebric multigrid PC - Mark Adams 2011
3  */
4 #include "petsc-private/matimpl.h"
5 #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6 #include <petsc-private/kspimpl.h>
7 
8 #if defined PETSC_GAMG_USE_LOG
9 PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10 #endif
11 
12 #if defined PETSC_USE_LOG
13 PetscLogEvent PC_GAMGGgraph_AGG;
14 PetscLogEvent PC_GAMGGgraph_GEO;
15 PetscLogEvent PC_GAMGCoarsen_AGG;
16 PetscLogEvent PC_GAMGCoarsen_GEO;
17 PetscLogEvent PC_GAMGProlongator_AGG;
18 PetscLogEvent PC_GAMGProlongator_GEO;
19 PetscLogEvent PC_GAMGOptprol_AGG;
20 PetscLogEvent PC_GAMGKKTProl_AGG;
21 #endif
22 
23 #define GAMG_MAXLEVELS 30
24 
25 /* #define GAMG_STAGES */
26 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
27 static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
28 #endif
29 
30 static PetscFList GAMGList = 0;
31 
32 /* ----------------------------------------------------------------------------- */
33 #undef __FUNCT__
34 #define __FUNCT__ "PCReset_GAMG"
35 PetscErrorCode PCReset_GAMG(PC pc)
36 {
37   PetscErrorCode  ierr;
38   PC_MG           *mg = (PC_MG*)pc->data;
39   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
40 
41   PetscFunctionBegin;
42   if( pc_gamg->data ) { /* this should not happen, cleaned up in SetUp */
43     ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr);
44   }
45   pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0;
46   PetscFunctionReturn(0);
47 }
48 
49 /* private 2x2 Mat Nest for Stokes */
50 typedef struct{
51   Mat A11,A21,A12,Amat;
52   IS  prim_is,constr_is;
53 }GAMGKKTMat;
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "GAMGKKTMatCreate"
57 static PetscErrorCode GAMGKKTMatCreate( Mat A, PetscBool iskkt, GAMGKKTMat *out )
58 {
59   PetscFunctionBegin;
60   out->Amat = A;
61   if( iskkt ){
62     PetscErrorCode   ierr;
63     IS       is_constraint, is_prime;
64     PetscInt nmin,nmax;
65 
66     ierr = MatGetOwnershipRange( A, &nmin, &nmax );   CHKERRQ(ierr);
67     ierr = MatFindZeroDiagonals( A, &is_constraint ); CHKERRQ(ierr);
68     ierr = ISComplement( is_constraint, nmin, nmax, &is_prime ); CHKERRQ(ierr);
69     out->prim_is = is_prime;
70     out->constr_is = is_constraint;
71 
72     ierr = MatGetSubMatrix( A, is_prime, is_prime,      MAT_INITIAL_MATRIX, &out->A11); CHKERRQ(ierr);
73     ierr = MatGetSubMatrix( A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12); CHKERRQ(ierr);
74     ierr = MatGetSubMatrix( A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21); CHKERRQ(ierr);
75 PetscPrintf(((PetscObject)A)->comm,"[%d]%s N=%d N_11=%d\n",0,__FUNCT__,A->rmap->N,out->A11->rmap->N);
76   }
77   else {
78     out->A11 = A;
79     out->A21 = PETSC_NULL;
80     out->A12 = PETSC_NULL;
81     out->prim_is = PETSC_NULL;
82     out->constr_is = PETSC_NULL;
83   }
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "GAMGKKTMatDestroy"
89 static PetscErrorCode GAMGKKTMatDestroy( GAMGKKTMat *mat )
90 {
91   PetscErrorCode   ierr;
92 
93   PetscFunctionBegin;
94   if( mat->A11 && mat->A11 != mat->Amat ) {
95     ierr = MatDestroy( &mat->A11 );  CHKERRQ(ierr);
96   }
97   ierr = MatDestroy( &mat->A21 );  CHKERRQ(ierr);
98   ierr = MatDestroy( &mat->A12 );  CHKERRQ(ierr);
99 
100   ierr = ISDestroy( &mat->prim_is );    CHKERRQ(ierr);
101   ierr = ISDestroy( &mat->constr_is );    CHKERRQ(ierr);
102 
103   PetscFunctionReturn(0);
104 }
105 
106 /* -------------------------------------------------------------------------- */
107 /*
108    createLevel: create coarse op with RAP.  repartition and/or reduce number
109      of active processors.
110 
111    Input Parameter:
112    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
113           'pc_gamg->data_sz' are changed via repartitioning/reduction.
114    . Amat_fine - matrix on this fine (k) level
115    . cr_bs - coarse block size
116    . isLast -
117    . stokes -
118    In/Output Parameter:
119    . a_P_inout - prolongation operator to the next level (k-->k-1)
120    . a_nactive_proc - number of active procs
121    Output Parameter:
122    . a_Amat_crs - coarse matrix that is created (k-1)
123 */
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "createLevel"
127 static PetscErrorCode createLevel( const PC pc,
128                                    const Mat Amat_fine,
129                                    const PetscInt cr_bs,
130                                    const PetscBool isLast,
131                                    const PetscBool stokes,
132                                    Mat *a_P_inout,
133                                    Mat *a_Amat_crs,
134                                    PetscMPIInt *a_nactive_proc
135                                    )
136 {
137   PetscErrorCode   ierr;
138   PC_MG           *mg = (PC_MG*)pc->data;
139   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
140   const PetscBool  repart = pc_gamg->repart;
141   const PetscInt   min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit;
142   Mat              Cmat,Pold=*a_P_inout;
143   MPI_Comm         wcomm = ((PetscObject)Amat_fine)->comm;
144   PetscMPIInt      mype,npe,new_npe,nactive=*a_nactive_proc;
145   PetscInt         ncrs_eq,ncrs_prim,f_bs;
146 
147   PetscFunctionBegin;
148   ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr);
149   ierr = MPI_Comm_size( wcomm, &npe );  CHKERRQ(ierr);
150   ierr = MatGetBlockSize( Amat_fine, &f_bs ); CHKERRQ(ierr);
151   /* RAP */
152   ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr);
153 
154   /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/
155   ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
156   assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0);
157   ierr = MatGetLocalSize( Cmat, &ncrs_eq, PETSC_NULL );  CHKERRQ(ierr);
158 
159   /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */
160   {
161     PetscInt ncrs_eq_glob;
162     ierr = MatGetSize( Cmat, &ncrs_eq_glob, PETSC_NULL );  CHKERRQ(ierr);
163     new_npe = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
164     if( new_npe == 0 || ncrs_eq_glob < coarse_max ) new_npe = 1;
165     else if ( new_npe >= nactive ) new_npe = nactive; /* no change, rare */
166     if( isLast ) new_npe = 1;
167   }
168 
169   if( !repart && new_npe==nactive ) {
170     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
171   }
172   else {
173     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
174     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
175     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
176     VecScatter      vecscat;
177     PetscScalar    *array;
178     Vec             src_crd, dest_crd;
179 
180     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
181 #if defined PETSC_GAMG_USE_LOG
182     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
183 #endif
184     /* make 'is_eq_newproc' */
185     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
186     if( repart && !stokes ) {
187       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
188       Mat             adj;
189 
190       if (pc_gamg->verbose>0) {
191         if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq);
192         else {
193           PetscInt n;
194           ierr = MPI_Allreduce( &ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm );CHKERRQ(ierr);
195           PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n);
196         }
197       }
198 
199       /* get 'adj' */
200       if( cr_bs == 1 ) {
201 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
202       }
203       else{
204 	/* make a scalar matrix to partition (no Stokes here) */
205 	Mat tMat;
206 	PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
207 	const PetscScalar *vals;
208 	const PetscInt *idx;
209 	PetscInt *d_nnz, *o_nnz, M, N;
210 	static PetscInt llev = 0;
211 
212 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
213 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
214         ierr = MatGetOwnershipRange( Cmat, &Istart_crs, &Iend_crs );    CHKERRQ(ierr);
215         ierr = MatGetSize( Cmat, &M, &N );CHKERRQ(ierr);
216 	for ( Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++ ) {
217 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
218 	  d_nnz[jj] = ncols/cr_bs;
219 	  o_nnz[jj] = ncols/cr_bs;
220 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
221 	  if( d_nnz[jj] > ncrs_prim ) d_nnz[jj] = ncrs_prim;
222 	  if( o_nnz[jj] > (M/cr_bs-ncrs_prim) ) o_nnz[jj] = M/cr_bs-ncrs_prim;
223 	}
224 
225 	ierr = MatCreate( wcomm, &tMat ); CHKERRQ(ierr);
226 	ierr = MatSetSizes( tMat, ncrs_prim, ncrs_prim,
227                             PETSC_DETERMINE, PETSC_DETERMINE );
228         CHKERRQ(ierr);
229         ierr = MatSetType(tMat,MATAIJ);   CHKERRQ(ierr);
230         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
231         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
232 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
233 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
234 
235 	for ( ii = Istart_crs; ii < Iend_crs; ii++ ) {
236 	  PetscInt dest_row = ii/cr_bs;
237 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
238 	  for( jj = 0 ; jj < ncols ; jj++ ){
239 	    PetscInt dest_col = idx[jj]/cr_bs;
240 	    PetscScalar v = 1.0;
241 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
242 	  }
243 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
244 	}
245 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
246 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247 
248 	if( llev++ == -1 ) {
249 	  PetscViewer viewer; char fname[32];
250 	  ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr);
251 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
252 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
253 	  ierr = PetscViewerDestroy( &viewer );
254 	}
255 
256 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
257 
258 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
259       } /* create 'adj' */
260 
261       { /* partition: get newproc_idx */
262 	char prefix[256];
263 	const char *pcpre;
264 	const PetscInt *is_idx;
265 	MatPartitioning  mpart;
266 	IS proc_is;
267 	PetscInt targetPE;
268 
269 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
270 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
271 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
272 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
273 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
274 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
275 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
276 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
277 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
278 
279 	/* collect IS info */
280 	ierr = PetscMalloc( ncrs_eq*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
281 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
282 	targetPE = 1; /* bring to "front" of machine */
283 	/*targetPE = npe/new_npe;*/ /* spread partitioning across machine */
284 	for( kk = jj = 0 ; kk < nloc_old ; kk++ ){
285 	  for( ii = 0 ; ii < cr_bs ; ii++, jj++ ){
286 	    newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
287 	  }
288 	}
289 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
290 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
291       }
292       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
293 
294       ierr = ISCreateGeneral( wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc );
295       CHKERRQ(ierr);
296       if( newproc_idx != 0 ) {
297 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
298       }
299     } /* repartitioning */
300     else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
301 
302       PetscInt rfactor,targetPE;
303       /* find factor */
304       if( new_npe == 1 ) rfactor = npe; /* easy */
305       else {
306 	PetscReal best_fact = 0.;
307 	jj = -1;
308 	for( kk = 1 ; kk <= npe ; kk++ ){
309 	  if( npe%kk==0 ) { /* a candidate */
310 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
311 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
312 	    if( fact > best_fact ) {
313 	      best_fact = fact; jj = kk;
314 	    }
315 	  }
316 	}
317 	if( jj != -1 ) rfactor = jj;
318 	else rfactor = 1; /* does this happen .. a prime */
319       }
320       new_npe = npe/rfactor;
321 
322       if( new_npe==nactive ) {
323 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
324 	ierr = PetscFree( counts );  CHKERRQ(ierr);
325 	if (pc_gamg->verbose>0){
326           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq);
327         }
328 	PetscFunctionReturn(0);
329       }
330 
331       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq);
332       targetPE = mype/rfactor;
333       ierr = ISCreateStride( wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc ); CHKERRQ(ierr);
334 
335       if( stokes ) {
336         ierr = ISCreateStride( wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim ); CHKERRQ(ierr);
337       }
338     } /* end simple 'is_eq_newproc' */
339 
340     /*
341      Create an index set from the is_eq_newproc index set to indicate the mapping TO
342      */
343     ierr = ISPartitioningToNumbering( is_eq_newproc, &is_eq_num ); CHKERRQ(ierr);
344     if( stokes ) {
345       ierr = ISPartitioningToNumbering( is_eq_newproc_prim, &is_eq_num_prim ); CHKERRQ(ierr);
346     }
347     else is_eq_num_prim = is_eq_num;
348     /*
349       Determine how many equations/vertices are assigned to each processor
350      */
351     ierr = ISPartitioningCount( is_eq_newproc, npe, counts ); CHKERRQ(ierr);
352     ncrs_eq_new = counts[mype];
353     ierr = ISDestroy( &is_eq_newproc );                       CHKERRQ(ierr);
354     if( stokes ) {
355       ierr = ISPartitioningCount( is_eq_newproc_prim, npe, counts ); CHKERRQ(ierr);
356       ierr = ISDestroy( &is_eq_newproc_prim );                       CHKERRQ(ierr);
357       ncrs_prim_new = counts[mype]/cr_bs; /* nodes */
358     }
359     else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
360 
361     ierr = PetscFree( counts );  CHKERRQ(ierr);
362 #if defined PETSC_GAMG_USE_LOG
363     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
364 #endif
365 
366     /* move data (for primal equations only) */
367     /* Create a vector to contain the newly ordered element information */
368     ierr = VecCreate( wcomm, &dest_crd );
369     ierr = VecSetSizes( dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE ); CHKERRQ(ierr);
370     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
371     /*
372      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
373      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
374      */
375     ierr = PetscMalloc( (ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
376     ierr = ISGetIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
377     for(ii=0,jj=0; ii<ncrs_prim ; ii++) {
378       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
379       for( kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
380     }
381     ierr = ISRestoreIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
382     ierr = ISCreateGeneral( wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat );
383     CHKERRQ(ierr);
384     ierr = PetscFree( tidx );  CHKERRQ(ierr);
385     /*
386      Create a vector to contain the original vertex information for each element
387      */
388     ierr = VecCreateSeq( PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd ); CHKERRQ(ierr);
389     for( jj=0; jj<ndata_cols ; jj++ ) {
390       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
391       for( ii=0 ; ii<ncrs_prim ; ii++) {
392 	for( kk=0; kk<ndata_rows ; kk++ ) {
393 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
394           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
395 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
396 	}
397       }
398     }
399     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
400     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
401     /*
402       Scatter the element vertex information (still in the original vertex ordering)
403       to the correct processor
404     */
405     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
406     CHKERRQ(ierr);
407     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
408     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
409     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
411     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
412     /*
413       Put the element vertex data into a new allocation of the gdata->ele
414     */
415     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
416     ierr = PetscMalloc( node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
417     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
418     strideNew = ncrs_prim_new*ndata_rows;
419     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
420     for( jj=0; jj<ndata_cols ; jj++ ) {
421       for( ii=0 ; ii<ncrs_prim_new ; ii++) {
422 	for( kk=0; kk<ndata_rows ; kk++ ) {
423 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
424 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
425 	}
426       }
427     }
428     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
429     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
430 
431     /* move A and P (columns) with new layout */
432 #if defined PETSC_GAMG_USE_LOG
433     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
434 #endif
435 
436     /*
437       Invert for MatGetSubMatrix
438     */
439     ierr = ISInvertPermutation( is_eq_num, ncrs_eq_new, &new_eq_indices ); CHKERRQ(ierr);
440     ierr = ISSort( new_eq_indices ); CHKERRQ(ierr); /* is this needed? */
441     ierr = ISSetBlockSize( new_eq_indices, cr_bs );   CHKERRQ(ierr);
442     if(is_eq_num != is_eq_num_prim) {
443       ierr = ISDestroy( &is_eq_num_prim ); CHKERRQ(ierr); /* could be same as 'is_eq_num' */
444     }
445     ierr = ISDestroy( &is_eq_num ); CHKERRQ(ierr);
446 #if defined PETSC_GAMG_USE_LOG
447     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
448     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
449 #endif
450     /* 'a_Amat_crs' output */
451     {
452       Mat mat;
453       ierr = MatGetSubMatrix( Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat );
454       CHKERRQ(ierr);
455       *a_Amat_crs = mat;
456 
457       if(!PETSC_TRUE){
458         PetscInt cbs, rbs;
459         ierr = MatGetBlockSizes( Cmat, &rbs, &cbs ); CHKERRQ(ierr);
460         PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
461         ierr = MatGetBlockSizes( mat, &rbs, &cbs ); CHKERRQ(ierr);
462         PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",mype,__FUNCT__,rbs,cbs,cr_bs);
463       }
464     }
465     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
466 
467 #if defined PETSC_GAMG_USE_LOG
468     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
469 #endif
470     /* prolongator */
471     {
472       IS findices;
473       PetscInt Istart,Iend;
474       Mat Pnew;
475       ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
476 #if defined PETSC_GAMG_USE_LOG
477       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
478 #endif
479       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
480       ierr = ISSetBlockSize(findices,f_bs);   CHKERRQ(ierr);
481       ierr = MatGetSubMatrix( Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew );
482       CHKERRQ(ierr);
483       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
484 
485       if(!PETSC_TRUE){
486         PetscInt cbs, rbs;
487         ierr = MatGetBlockSizes( Pold, &rbs, &cbs ); CHKERRQ(ierr);
488         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
489         ierr = MatGetBlockSizes( Pnew, &rbs, &cbs ); CHKERRQ(ierr);
490         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
491       }
492 #if defined PETSC_GAMG_USE_LOG
493       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
494 #endif
495       ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
496 
497       /* output - repartitioned */
498       *a_P_inout = Pnew;
499     }
500     ierr = ISDestroy( &new_eq_indices ); CHKERRQ(ierr);
501 
502     *a_nactive_proc = new_npe; /* output */
503   }
504 
505   /* outout matrix data */
506   if( !PETSC_TRUE ) {
507     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
508     if(llev==0) {
509       sprintf(fname,"Cmat_%d.m",llev++);
510       PetscViewerASCIIOpen(wcomm,fname,&viewer);
511       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
512       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
513       ierr = PetscViewerDestroy( &viewer );
514     }
515     sprintf(fname,"Cmat_%d.m",llev++);
516     PetscViewerASCIIOpen(wcomm,fname,&viewer);
517     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
518     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
519     ierr = PetscViewerDestroy( &viewer );
520   }
521 
522   PetscFunctionReturn(0);
523 }
524 
525 /* -------------------------------------------------------------------------- */
526 /*
527    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
528                     by setting data structures and options.
529 
530    Input Parameter:
531 .  pc - the preconditioner context
532 
533    Application Interface Routine: PCSetUp()
534 
535    Notes:
536    The interface routine PCSetUp() is not usually called directly by
537    the user, but instead is called by PCApply() if necessary.
538 */
539 #undef __FUNCT__
540 #define __FUNCT__ "PCSetUp_GAMG"
541 PetscErrorCode PCSetUp_GAMG( PC pc )
542 {
543   PetscErrorCode  ierr;
544   PC_MG           *mg = (PC_MG*)pc->data;
545   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
546   Mat              Pmat = pc->pmat;
547   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
548   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
549   PetscMPIInt      mype,npe,nactivepe;
550   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
551   PetscReal        emaxs[GAMG_MAXLEVELS];
552   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS];
553   PetscInt         level_bs[GAMG_MAXLEVELS];
554   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
555   PetscLogDouble   nnz0=0.,nnztot=0.;
556   MatInfo          info;
557   PetscBool        stokes = PETSC_FALSE;
558 
559   PetscFunctionBegin;
560   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
561   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
562   if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled);
563   if( pc_gamg->setup_count++ > 0 ) {
564     PC_MG_Levels **mglevels = mg->levels;
565     /* just do Galerkin grids */
566     Mat B,dA,dB;
567     assert(pc->setupcalled);
568 
569     if( pc_gamg->Nlevels > 1 ) {
570       /* currently only handle case where mat and pmat are the same on coarser levels */
571       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
572       /* (re)set to get dirty flag */
573       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
574 
575       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
576         /* the first time through the matrix structure has changed from repartitioning */
577         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
578           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
579           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
580           mglevels[level]->A = B;
581         }
582         else {
583           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
584           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
585         }
586         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
587         dB = B;
588       }
589     }
590 
591     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
592 
593     /* PCSetUp_MG seems to insists on setting this to GMRES */
594     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
595 
596     PetscFunctionReturn(0);
597   }
598   assert(pc->setupcalled == 0);
599 
600   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
601 
602   ierr = GAMGKKTMatCreate( Pmat, stokes, &kktMatsArr[0] ); CHKERRQ(ierr);
603 
604   if( pc_gamg->data == 0 ) {
605     if( !pc_gamg->createdefaultdata ){
606       SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
607     }
608     if( stokes ) {
609       SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
610     }
611     ierr = pc_gamg->createdefaultdata( pc, kktMatsArr[0].A11 ); CHKERRQ(ierr);
612   }
613 
614   /* get basic dims */
615   if( stokes ) {
616     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
617   }
618   else {
619     ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
620   }
621 
622   ierr = MatGetSize( Pmat, &M, &qq );CHKERRQ(ierr);
623   if (pc_gamg->verbose) {
624     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
625     else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
626     CHKERRQ(ierr);
627     nnz0 = info.nz_used;
628     nnztot = info.nz_used;
629     PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
630                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
631                 (int)(nnz0/(PetscReal)M),npe);
632   }
633 
634   /* Get A_i and R_i */
635   for ( level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
636         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
637         level++ ){
638     level1 = level + 1;
639 #if defined PETSC_GAMG_USE_LOG
640     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
641 #if (defined GAMG_STAGES)
642     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
643 #endif
644 #endif
645     /* deal with Stokes, get sub matrices */
646     if( level > 0 ) {
647       ierr = GAMGKKTMatCreate( Aarr[level], stokes, &kktMatsArr[level] ); CHKERRQ(ierr);
648     }
649     { /* construct prolongator */
650       Mat Gmat;
651       PetscCoarsenData *agg_lists;
652       Mat Prol11,Prol22;
653 
654       level_bs[level] = bs;
655       ierr = pc_gamg->graph( pc,kktMatsArr[level].A11, &Gmat ); CHKERRQ(ierr);
656       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
657       ierr = pc_gamg->prolongator( pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11 ); CHKERRQ(ierr);
658 
659       /* could have failed to create new level */
660       if( Prol11 ){
661         /* get new block size of coarse matrices */
662         ierr = MatGetBlockSizes( Prol11, PETSC_NULL, &bs ); CHKERRQ(ierr);
663 
664         if( stokes ) {
665           if(!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
666           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
667           ierr = pc_gamg->formkktprol( pc, Prol11, kktMatsArr[level].A21, &Prol22 ); CHKERRQ(ierr);
668         }
669 
670         if( pc_gamg->optprol ){
671           /* smooth */
672           ierr = pc_gamg->optprol( pc, kktMatsArr[level].A11, &Prol11 ); CHKERRQ(ierr);
673         }
674 
675         /* remove rows of singleton rows (BCs) */
676 
677         if( stokes ) {
678           IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is};
679           Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 };
680           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1] ); CHKERRQ(ierr);
681         }
682         else {
683           Parr[level1] = Prol11;
684         }
685       }
686       else Parr[level1] = PETSC_NULL;
687 
688       if ( pc_gamg->use_aggs_in_gasm ) {
689         ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level] );  CHKERRQ(ierr);
690       }
691 
692       ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] );  CHKERRQ(ierr);
693 
694       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
695       ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
696     } /* construct prolongator scope */
697 #if defined PETSC_GAMG_USE_LOG
698     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
699 #endif
700     /* cache eigen estimate */
701     if( pc_gamg->emax_id != -1 ){
702       PetscBool flag;
703       ierr = PetscObjectComposedDataGetReal( (PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag );
704       CHKERRQ( ierr );
705       if( !flag ) emaxs[level] = -1.;
706     }
707     else emaxs[level] = -1.;
708     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
709     if( !Parr[level1] ) {
710       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
711       break;
712     }
713 #if defined PETSC_GAMG_USE_LOG
714     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
715 #endif
716 
717     ierr = createLevel( pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
718                         stokes, &Parr[level1], &Aarr[level1], &nactivepe );
719     CHKERRQ(ierr);
720 
721 #if defined PETSC_GAMG_USE_LOG
722     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
723 #endif
724     ierr = MatGetSize( Aarr[level1], &M, &qq );CHKERRQ(ierr);
725 
726     if (pc_gamg->verbose > 0){
727       PetscInt NN = M;
728       if(pc_gamg->verbose==1) {
729         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); CHKERRQ(ierr);
730         ierr = MatGetLocalSize( Aarr[level1], &NN, &qq );
731       }
732       else ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info );
733 
734       CHKERRQ(ierr);
735       nnztot += info.nz_used;
736       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
737                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
738                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
739       CHKERRQ(ierr);
740     }
741 
742     /* stop if one node -- could pull back for singular problems */
743     if( M/pc_gamg->data_cell_cols < 2 ) {
744       level++;
745       break;
746     }
747 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
748     ierr = PetscLogStagePop(); CHKERRQ( ierr );
749 #endif
750   } /* levels */
751 
752   if( pc_gamg->data ) {
753     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
754     pc_gamg->data = PETSC_NULL;
755   }
756 
757   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
758   pc_gamg->Nlevels = level + 1;
759   fine_level = level;
760   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
761 
762   /* simple setup */
763   if( !PETSC_TRUE ){
764     PC_MG_Levels **mglevels = mg->levels;
765     for (lidx=0,level=pc_gamg->Nlevels-1;
766          lidx<fine_level;
767          lidx++, level--){
768       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
769       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr);
770       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
771       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
772     }
773     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
774 
775     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
776   }
777   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
778     /* set default smoothers & set operators */
779     for ( lidx = 1, level = pc_gamg->Nlevels-2;
780           lidx <= fine_level;
781           lidx++, level--) {
782       KSP smoother;
783       PC subpc;
784 
785       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
786       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
787 
788       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
789       /* set ops */
790       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
791       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
792 
793       /* create field split PC, get subsmoother */
794       if( stokes ) {
795         KSP *ksps;
796         PetscInt nn;
797         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);   CHKERRQ(ierr);
798         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is); CHKERRQ(ierr);
799         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
800         smoother = ksps[0];
801         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
802         ierr = PetscFree( ksps ); CHKERRQ(ierr);
803       }
804       ierr = GAMGKKTMatDestroy( &kktMatsArr[level] ); CHKERRQ(ierr);
805 
806       /* set defaults */
807       ierr = KSPSetType( smoother, KSPCHEBYSHEV );CHKERRQ(ierr);
808 
809       /* override defaults and command line args (!) */
810       if ( pc_gamg->use_aggs_in_gasm ) {
811         PetscInt sz;
812         IS *is;
813 
814         sz = nASMBlocksArr[level];
815         is = ASMLocalIDsArr[level];
816         ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr);
817         if(sz==0){
818           IS is;
819           PetscInt my0,kk;
820           ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr);
821           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr);
822           ierr = PCGASMSetSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr);
823           ierr = ISDestroy( &is ); CHKERRQ(ierr);
824         }
825         else {
826           PetscInt kk;
827           ierr = PCGASMSetSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr);
828           for(kk=0;kk<sz;kk++){
829             ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr);
830           }
831           ierr = PetscFree( is ); CHKERRQ(ierr);
832         }
833         ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr);
834 
835         ASMLocalIDsArr[level] = PETSC_NULL;
836         nASMBlocksArr[level] = 0;
837         ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr);
838       }
839       else {
840         ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
841       }
842     }
843     {
844       /* coarse grid */
845       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
846       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
847       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
848       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
849       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
850       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
851       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
852       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
853       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
854       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
855       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
856     }
857 
858     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
859     ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr);
860     ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr);
861     ierr = PetscOptionsEnd();  CHKERRQ(ierr);
862     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
863 
864     /* create cheby smoothers */
865     for ( lidx = 1, level = pc_gamg->Nlevels-2;
866           lidx <= fine_level;
867           lidx++, level--) {
868       KSP smoother;
869       PetscBool flag;
870       PC subpc;
871 
872       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
873       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
874 
875       /* create field split PC, get subsmoother */
876       if( stokes ) {
877         KSP *ksps;
878         PetscInt nn;
879         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
880         smoother = ksps[0];
881         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
882         ierr = PetscFree( ksps );  CHKERRQ(ierr);
883       }
884 
885       /* do my own cheby */
886       ierr = PetscObjectTypeCompare( (PetscObject)smoother, KSPCHEBYSHEV, &flag ); CHKERRQ(ierr);
887       if( flag ) {
888         PetscReal emax, emin;
889         ierr = PetscObjectTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr);
890         if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
891         else{ /* eigen estimate 'emax' */
892           KSP eksp; Mat Lmat = Aarr[level];
893           Vec bb, xx;
894 
895           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
896           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
897           {
898             PetscRandom    rctx;
899             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
900             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
901             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
902             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
903           }
904 
905           if( removedEqs[level] ) {
906             /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */
907             PetscScalar *zeros;
908             PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level];
909             const PetscInt *idx;
910             ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr);
911             ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr);
912             for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.;
913             ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr);
914             ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr);
915             for(ii=0;ii<sz;ii++) {
916               for(jj=0;jj<bs;jj++) {
917                 idx_bs[ii] = bs*idx[ii]+jj;
918               }
919             }
920             ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr);
921             if( sz > 0 ) {
922               ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES );  CHKERRQ(ierr);
923             }
924             ierr = PetscFree( idx_bs );  CHKERRQ(ierr);
925             ierr = PetscFree( zeros );  CHKERRQ(ierr);
926             ierr = VecAssemblyBegin(bb); CHKERRQ(ierr);
927             ierr = VecAssemblyEnd(bb); CHKERRQ(ierr);
928           }
929           ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr);
930           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
931           CHKERRQ(ierr);
932           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
933           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
934           ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_");         CHKERRQ(ierr);
935           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
936 
937           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
938           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
939           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
940 
941           /* set PC type to be same as smoother */
942           ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr );
943 
944           /* solve - keep stuff out of logging */
945           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
946           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
947           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
948           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
949           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
950 
951           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
952 
953           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
954           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
955           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
956 
957           if( pc_gamg->verbose > 0 ) {
958             PetscInt N1, tt;
959             ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
960             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
961           }
962         }
963         {
964           PetscInt N1, N0;
965           ierr = MatGetSize( Aarr[level], &N1, PETSC_NULL );         CHKERRQ(ierr);
966           ierr = MatGetSize( Aarr[level+1], &N0, PETSC_NULL );       CHKERRQ(ierr);
967           /* heuristic - is this crap? */
968           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
969 	  emin = emax * pc_gamg->eigtarget[0];
970           emax *= pc_gamg->eigtarget[1];
971         }
972         ierr = KSPChebyshevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr);
973       } /* setup checby flag */
974 
975       if( removedEqs[level] ) {
976         ierr = ISDestroy( &removedEqs[level] );                    CHKERRQ(ierr);
977       }
978     } /* non-coarse levels */
979 
980     /* clean up */
981     for(level=1;level<pc_gamg->Nlevels;level++){
982       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
983       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
984     }
985 
986     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
987 
988     if( PETSC_FALSE ){
989       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
990       ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
991       ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
992     }
993   }
994   else {
995     KSP smoother;
996     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
997     ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr);
998     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
999     ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr);
1000     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
1001   }
1002 
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 /* ------------------------------------------------------------------------- */
1007 /*
1008  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
1009    that was created with PCCreate_GAMG().
1010 
1011    Input Parameter:
1012 .  pc - the preconditioner context
1013 
1014    Application Interface Routine: PCDestroy()
1015 */
1016 #undef __FUNCT__
1017 #define __FUNCT__ "PCDestroy_GAMG"
1018 PetscErrorCode PCDestroy_GAMG( PC pc )
1019 {
1020   PetscErrorCode  ierr;
1021   PC_MG           *mg = (PC_MG*)pc->data;
1022   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
1023 
1024   PetscFunctionBegin;
1025   ierr = PCReset_GAMG( pc );CHKERRQ(ierr);
1026   ierr = PetscFree( pc_gamg );CHKERRQ(ierr);
1027   ierr = PCDestroy_MG( pc );CHKERRQ(ierr);
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 
1032 #undef __FUNCT__
1033 #define __FUNCT__ "PCGAMGSetProcEqLim"
1034 /*@
1035    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1036    processor reduction.
1037 
1038    Not Collective on PC
1039 
1040    Input Parameters:
1041 .  pc - the preconditioner context
1042 
1043 
1044    Options Database Key:
1045 .  -pc_gamg_process_eq_limit
1046 
1047    Level: intermediate
1048 
1049    Concepts: Unstructured multrigrid preconditioner
1050 
1051 .seealso: ()
1052 @*/
1053 PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1054 {
1055   PetscErrorCode ierr;
1056 
1057   PetscFunctionBegin;
1058   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1059   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 EXTERN_C_BEGIN
1064 #undef __FUNCT__
1065 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1066 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1067 {
1068   PC_MG           *mg = (PC_MG*)pc->data;
1069   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1070 
1071   PetscFunctionBegin;
1072   if(n>0) pc_gamg->min_eq_proc = n;
1073   PetscFunctionReturn(0);
1074 }
1075 EXTERN_C_END
1076 
1077 #undef __FUNCT__
1078 #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1079 /*@
1080    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1081 
1082  Collective on PC
1083 
1084    Input Parameters:
1085 .  pc - the preconditioner context
1086 
1087 
1088    Options Database Key:
1089 .  -pc_gamg_coarse_eq_limit
1090 
1091    Level: intermediate
1092 
1093    Concepts: Unstructured multrigrid preconditioner
1094 
1095 .seealso: ()
1096  @*/
1097 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1098 {
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1103   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 EXTERN_C_BEGIN
1108 #undef __FUNCT__
1109 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1110 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1111 {
1112   PC_MG           *mg = (PC_MG*)pc->data;
1113   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1114 
1115   PetscFunctionBegin;
1116   if(n>0) pc_gamg->coarse_eq_limit = n;
1117   PetscFunctionReturn(0);
1118 }
1119 EXTERN_C_END
1120 
1121 #undef __FUNCT__
1122 #define __FUNCT__ "PCGAMGSetRepartitioning"
1123 /*@
1124    PCGAMGSetRepartitioning - Repartition the coarse grids
1125 
1126    Collective on PC
1127 
1128    Input Parameters:
1129 .  pc - the preconditioner context
1130 
1131 
1132    Options Database Key:
1133 .  -pc_gamg_repartition
1134 
1135    Level: intermediate
1136 
1137    Concepts: Unstructured multrigrid preconditioner
1138 
1139 .seealso: ()
1140 @*/
1141 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1142 {
1143   PetscErrorCode ierr;
1144 
1145   PetscFunctionBegin;
1146   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1147   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 EXTERN_C_BEGIN
1152 #undef __FUNCT__
1153 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
1154 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1155 {
1156   PC_MG           *mg = (PC_MG*)pc->data;
1157   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1158 
1159   PetscFunctionBegin;
1160   pc_gamg->repart = n;
1161   PetscFunctionReturn(0);
1162 }
1163 EXTERN_C_END
1164 
1165 #undef __FUNCT__
1166 #define __FUNCT__ "PCGAMGSetUseASMAggs"
1167 /*@
1168    PCGAMGSetUseASMAggs -
1169 
1170    Collective on PC
1171 
1172    Input Parameters:
1173 .  pc - the preconditioner context
1174 
1175 
1176    Options Database Key:
1177 .  -pc_gamg_use_agg_gasm
1178 
1179    Level: intermediate
1180 
1181    Concepts: Unstructured multrigrid preconditioner
1182 
1183 .seealso: ()
1184 @*/
1185 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1186 {
1187   PetscErrorCode ierr;
1188 
1189   PetscFunctionBegin;
1190   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1191   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1192   PetscFunctionReturn(0);
1193 }
1194 
1195 EXTERN_C_BEGIN
1196 #undef __FUNCT__
1197 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1198 PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1199 {
1200   PC_MG           *mg = (PC_MG*)pc->data;
1201   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1202 
1203   PetscFunctionBegin;
1204   pc_gamg->use_aggs_in_gasm = n;
1205   PetscFunctionReturn(0);
1206 }
1207 EXTERN_C_END
1208 
1209 #undef __FUNCT__
1210 #define __FUNCT__ "PCGAMGSetNlevels"
1211 /*@
1212    PCGAMGSetNlevels -
1213 
1214    Not collective on PC
1215 
1216    Input Parameters:
1217 .  pc - the preconditioner context
1218 
1219 
1220    Options Database Key:
1221 .  -pc_mg_levels
1222 
1223    Level: intermediate
1224 
1225    Concepts: Unstructured multrigrid preconditioner
1226 
1227 .seealso: ()
1228 @*/
1229 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1230 {
1231   PetscErrorCode ierr;
1232 
1233   PetscFunctionBegin;
1234   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1235   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 EXTERN_C_BEGIN
1240 #undef __FUNCT__
1241 #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
1242 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1243 {
1244   PC_MG           *mg = (PC_MG*)pc->data;
1245   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1246 
1247   PetscFunctionBegin;
1248   pc_gamg->Nlevels = n;
1249   PetscFunctionReturn(0);
1250 }
1251 EXTERN_C_END
1252 
1253 #undef __FUNCT__
1254 #define __FUNCT__ "PCGAMGSetThreshold"
1255 /*@
1256    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1257 
1258    Not collective on PC
1259 
1260    Input Parameters:
1261 .  pc - the preconditioner context
1262 
1263 
1264    Options Database Key:
1265 .  -pc_gamg_threshold
1266 
1267    Level: intermediate
1268 
1269    Concepts: Unstructured multrigrid preconditioner
1270 
1271 .seealso: ()
1272 @*/
1273 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
1274 {
1275   PetscErrorCode ierr;
1276 
1277   PetscFunctionBegin;
1278   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1279   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
1280   PetscFunctionReturn(0);
1281 }
1282 
1283 EXTERN_C_BEGIN
1284 #undef __FUNCT__
1285 #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
1286 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
1287 {
1288   PC_MG           *mg = (PC_MG*)pc->data;
1289   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1290 
1291   PetscFunctionBegin;
1292   pc_gamg->threshold = n;
1293   PetscFunctionReturn(0);
1294 }
1295 EXTERN_C_END
1296 
1297 #undef __FUNCT__
1298 #define __FUNCT__ "PCGAMGSetType"
1299 /*@
1300    PCGAMGSetType - Set solution method - calls sub create method
1301 
1302    Collective on PC
1303 
1304    Input Parameters:
1305 .  pc - the preconditioner context
1306 
1307 
1308    Options Database Key:
1309 .  -pc_gamg_type
1310 
1311    Level: intermediate
1312 
1313    Concepts: Unstructured multrigrid preconditioner
1314 
1315 .seealso: ()
1316 @*/
1317 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type )
1318 {
1319   PetscErrorCode ierr;
1320 
1321   PetscFunctionBegin;
1322   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1323   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type));
1324   CHKERRQ(ierr);
1325   PetscFunctionReturn(0);
1326 }
1327 
1328 EXTERN_C_BEGIN
1329 #undef __FUNCT__
1330 #define __FUNCT__ "PCGAMGSetType_GAMG"
1331 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type )
1332 {
1333   PetscErrorCode ierr,(*r)(PC);
1334 
1335   PetscFunctionBegin;
1336   ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r);
1337   CHKERRQ(ierr);
1338 
1339   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1340 
1341   /* call sub create method */
1342   ierr = (*r)(pc); CHKERRQ(ierr);
1343 
1344   PetscFunctionReturn(0);
1345 }
1346 EXTERN_C_END
1347 
1348 #undef __FUNCT__
1349 #define __FUNCT__ "PCSetFromOptions_GAMG"
1350 PetscErrorCode PCSetFromOptions_GAMG( PC pc )
1351 {
1352   PetscErrorCode  ierr;
1353   PC_MG           *mg = (PC_MG*)pc->data;
1354   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1355   PetscBool        flag;
1356   PetscInt         two = 2;
1357   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1358 
1359   PetscFunctionBegin;
1360   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1361   {
1362     /* -pc_gamg_verbose */
1363     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
1364                            "none", pc_gamg->verbose,
1365                            &pc_gamg->verbose, PETSC_NULL );
1366     CHKERRQ(ierr);
1367 
1368     /* -pc_gamg_repartition */
1369     ierr = PetscOptionsBool("-pc_gamg_repartition",
1370                             "Repartion coarse grids (false)",
1371                             "PCGAMGRepartitioning",
1372                             pc_gamg->repart,
1373                             &pc_gamg->repart,
1374                             &flag);
1375     CHKERRQ(ierr);
1376 
1377     /* -pc_gamg_use_agg_gasm */
1378     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1379                             "Use aggregation agragates for GASM smoother (false)",
1380                             "PCGAMGUseASMAggs",
1381                             pc_gamg->use_aggs_in_gasm,
1382                             &pc_gamg->use_aggs_in_gasm,
1383                             &flag);
1384     CHKERRQ(ierr);
1385 
1386     /* -pc_gamg_process_eq_limit */
1387     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1388                            "Limit (goal) on number of equations per process on coarse grids",
1389                            "PCGAMGSetProcEqLim",
1390                            pc_gamg->min_eq_proc,
1391                            &pc_gamg->min_eq_proc,
1392                            &flag );
1393     CHKERRQ(ierr);
1394 
1395     /* -pc_gamg_coarse_eq_limit */
1396     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1397                            "Limit on number of equations for the coarse grid",
1398                            "PCGAMGSetCoarseEqLim",
1399                            pc_gamg->coarse_eq_limit,
1400                            &pc_gamg->coarse_eq_limit,
1401                            &flag );
1402     CHKERRQ(ierr);
1403 
1404     /* -pc_gamg_threshold */
1405     ierr = PetscOptionsReal("-pc_gamg_threshold",
1406                             "Relative threshold to use for dropping edges in aggregation graph",
1407                             "PCGAMGSetThreshold",
1408                             pc_gamg->threshold,
1409                             &pc_gamg->threshold,
1410                             &flag );
1411     CHKERRQ(ierr);
1412     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
1413 
1414     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
1415 
1416     ierr = PetscOptionsInt("-pc_mg_levels",
1417                            "Set number of MG levels",
1418                            "PCGAMGSetNlevels",
1419                            pc_gamg->Nlevels,
1420                            &pc_gamg->Nlevels,
1421                            &flag );
1422   }
1423   ierr = PetscOptionsTail();CHKERRQ(ierr);
1424 
1425   PetscFunctionReturn(0);
1426 }
1427 
1428 /* -------------------------------------------------------------------------- */
1429 /*MC
1430      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
1431        - This is the entry point to GAMG, registered in pcregis.c
1432 
1433    Options Database Keys:
1434    Multigrid options(inherited)
1435 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1436 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1437 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1438 -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1439 
1440   Level: intermediate
1441 
1442   Concepts: multigrid
1443 
1444 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1445            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1446            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1447            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1448            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1449 M*/
1450 EXTERN_C_BEGIN
1451 #undef __FUNCT__
1452 #define __FUNCT__ "PCCreate_GAMG"
1453 PetscErrorCode  PCCreate_GAMG( PC pc )
1454 {
1455   PetscErrorCode  ierr;
1456   PC_GAMG         *pc_gamg;
1457   PC_MG           *mg;
1458 #if defined PETSC_GAMG_USE_LOG
1459   static long count = 0;
1460 #endif
1461 
1462   PetscFunctionBegin;
1463 
1464   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1465   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
1466   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
1467 
1468   /* create a supporting struct and attach it to pc */
1469   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
1470   mg = (PC_MG*)pc->data;
1471   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
1472   mg->innerctx = pc_gamg;
1473 
1474   pc_gamg->setup_count = 0;
1475   /* these should be in subctx but repartitioning needs simple arrays */
1476   pc_gamg->data_sz = 0;
1477   pc_gamg->data = 0;
1478 
1479   /* register AMG type */
1480   if( !GAMGList ){
1481     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1482     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
1483   }
1484 
1485   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1486   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1487   pc->ops->setup          = PCSetUp_GAMG;
1488   pc->ops->reset          = PCReset_GAMG;
1489   pc->ops->destroy        = PCDestroy_GAMG;
1490 
1491   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1492 					    "PCGAMGSetProcEqLim_C",
1493 					    "PCGAMGSetProcEqLim_GAMG",
1494 					    PCGAMGSetProcEqLim_GAMG);
1495   CHKERRQ(ierr);
1496 
1497   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1498 					    "PCGAMGSetCoarseEqLim_C",
1499 					    "PCGAMGSetCoarseEqLim_GAMG",
1500 					    PCGAMGSetCoarseEqLim_GAMG);
1501   CHKERRQ(ierr);
1502 
1503   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1504 					    "PCGAMGSetRepartitioning_C",
1505 					    "PCGAMGSetRepartitioning_GAMG",
1506 					    PCGAMGSetRepartitioning_GAMG);
1507   CHKERRQ(ierr);
1508 
1509   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1510 					    "PCGAMGSetUseASMAggs_C",
1511 					    "PCGAMGSetUseASMAggs_GAMG",
1512 					    PCGAMGSetUseASMAggs_GAMG);
1513   CHKERRQ(ierr);
1514 
1515   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1516 					    "PCGAMGSetThreshold_C",
1517 					    "PCGAMGSetThreshold_GAMG",
1518 					    PCGAMGSetThreshold_GAMG);
1519   CHKERRQ(ierr);
1520 
1521   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1522 					    "PCGAMGSetType_C",
1523 					    "PCGAMGSetType_GAMG",
1524 					    PCGAMGSetType_GAMG);
1525   CHKERRQ(ierr);
1526 
1527   pc_gamg->repart = PETSC_FALSE;
1528   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1529   pc_gamg->min_eq_proc = 100;
1530   pc_gamg->coarse_eq_limit = 800;
1531   pc_gamg->threshold = 0.001;
1532   pc_gamg->Nlevels = GAMG_MAXLEVELS;
1533   pc_gamg->verbose = 0;
1534   pc_gamg->emax_id = -1;
1535   pc_gamg->eigtarget[0] = 0.05;
1536   pc_gamg->eigtarget[1] = 1.05;
1537 
1538   /* private events */
1539 #if defined PETSC_GAMG_USE_LOG
1540   if( count++ == 0 ) {
1541     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1542     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1543     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1544     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1545     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1546     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1547     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1548     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1549     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1550     PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1551     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1552     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1553     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1554     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1555     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1556     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1557     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1558 
1559     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1560     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1561     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1562     /* create timer stages */
1563 #if defined GAMG_STAGES
1564     {
1565       char str[32];
1566       sprintf(str,"MG Level %d (finest)",0);
1567       PetscLogStageRegister(str, &gamg_stages[0]);
1568       PetscInt lidx;
1569       for (lidx=1;lidx<9;lidx++){
1570 	sprintf(str,"MG Level %d",lidx);
1571 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1572       }
1573     }
1574 #endif
1575   }
1576 #endif
1577   /* general events */
1578 #if defined PETSC_USE_LOG
1579   PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);
1580   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
1581   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1582   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1583   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1584   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1585   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
1586   PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);
1587 #endif
1588 
1589   /* instantiate derived type */
1590   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1591   {
1592     char tname[256] = GAMGAGG;
1593     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
1594                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
1595     CHKERRQ(ierr);
1596     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
1597   }
1598   ierr = PetscOptionsTail();   CHKERRQ(ierr);
1599 
1600   PetscFunctionReturn(0);
1601 }
1602 EXTERN_C_END
1603