xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 596f5af7c3d0c0d6c3f65e9efcda2bd4733ba063)
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,ncrs_eq_ave;
162     ierr = MatGetSize( Cmat, &ncrs_eq_glob, PETSC_NULL );  CHKERRQ(ierr);
163     ncrs_eq_ave = ncrs_eq_glob/npe;
164     new_npe = (PetscMPIInt)((float)ncrs_eq_ave/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
165     if( new_npe == 0 || ncrs_eq_ave < coarse_max ) new_npe = 1;
166     else if ( new_npe >= nactive ) new_npe = nactive; /* no change, rare */
167     if( isLast ) new_npe = 1;
168   }
169 
170   if( !repart && new_npe==nactive ) {
171     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
172     }
173   else {
174     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
175     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old;
176     IS              is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices;
177     VecScatter      vecscat;
178     PetscScalar    *array;
179     Vec             src_crd, dest_crd;
180 
181     nloc_old = ncrs_eq/cr_bs; assert(ncrs_eq%cr_bs==0);
182 #if defined PETSC_GAMG_USE_LOG
183     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
184 #endif
185     /* make 'is_eq_newproc' */
186     ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr);
187     if( repart && !stokes ) {
188       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
189       Mat             adj;
190 
191       if (pc_gamg->verbose>0) {
192         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);
193         else {
194           PetscInt n;
195           ierr = MPI_Allreduce( &ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm );CHKERRQ(ierr);
196           PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n);
197         }
198       }
199 
200       /* get 'adj' */
201       if( cr_bs == 1 ) {
202 	ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
203       }
204       else{
205 	/* make a scalar matrix to partition (no Stokes here) */
206 	Mat tMat;
207 	PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
208 	const PetscScalar *vals;
209 	const PetscInt *idx;
210 	PetscInt *d_nnz, *o_nnz, M, N;
211 	static PetscInt llev = 0;
212 
213 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr);
214 	ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr);
215         ierr = MatGetOwnershipRange( Cmat, &Istart_crs, &Iend_crs );    CHKERRQ(ierr);
216         ierr = MatGetSize( Cmat, &M, &N );CHKERRQ(ierr);
217 	for ( Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cr_bs, jj++ ) {
218 	  ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
219 	  d_nnz[jj] = ncols/cr_bs;
220 	  o_nnz[jj] = ncols/cr_bs;
221 	  ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr);
222 	  if( d_nnz[jj] > ncrs_prim ) d_nnz[jj] = ncrs_prim;
223 	  if( o_nnz[jj] > (M/cr_bs-ncrs_prim) ) o_nnz[jj] = M/cr_bs-ncrs_prim;
224 	}
225 
226 	ierr = MatCreate( wcomm, &tMat ); CHKERRQ(ierr);
227 	ierr = MatSetSizes( tMat, ncrs_prim, ncrs_prim,
228                             PETSC_DETERMINE, PETSC_DETERMINE );
229         CHKERRQ(ierr);
230         ierr = MatSetType(tMat,MATAIJ);   CHKERRQ(ierr);
231         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
232         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
233 	ierr = PetscFree( d_nnz ); CHKERRQ(ierr);
234 	ierr = PetscFree( o_nnz ); CHKERRQ(ierr);
235 
236 	for ( ii = Istart_crs; ii < Iend_crs; ii++ ) {
237 	  PetscInt dest_row = ii/cr_bs;
238 	  ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
239 	  for( jj = 0 ; jj < ncols ; jj++ ){
240 	    PetscInt dest_col = idx[jj]/cr_bs;
241 	    PetscScalar v = 1.0;
242 	    ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr);
243 	  }
244 	  ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr);
245 	}
246 	ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247 	ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248 
249 	if( llev++ == -1 ) {
250 	  PetscViewer viewer; char fname[32];
251 	  ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr);
252 	  PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer);
253 	  ierr = MatView( tMat, viewer ); CHKERRQ(ierr);
254 	  ierr = PetscViewerDestroy( &viewer );
255 	}
256 
257 	ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj );   CHKERRQ(ierr);
258 
259 	ierr = MatDestroy( &tMat );  CHKERRQ(ierr);
260       } /* create 'adj' */
261 
262       { /* partition: get newproc_idx */
263 	char prefix[256];
264 	const char *pcpre;
265 	const PetscInt *is_idx;
266 	MatPartitioning  mpart;
267 	IS proc_is;
268 	PetscInt targetPE;
269 
270 	ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr);
271 	ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr);
272 	ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr);
273 	ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr);
274 	ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
275 	ierr = MatPartitioningSetFromOptions( mpart );    CHKERRQ(ierr);
276 	ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr);
277 	ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr);
278 	ierr = MatPartitioningDestroy( &mpart );          CHKERRQ(ierr);
279 
280 	/* collect IS info */
281 	ierr = PetscMalloc( ncrs_eq*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr);
282 	ierr = ISGetIndices( proc_is, &is_idx );        CHKERRQ(ierr);
283 	targetPE = 1; /* bring to "front" of machine */
284 	/*targetPE = npe/new_npe;*/ /* spread partitioning across machine */
285 	for( kk = jj = 0 ; kk < nloc_old ; kk++ ){
286 	  for( ii = 0 ; ii < cr_bs ; ii++, jj++ ){
287 	    newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
288 	  }
289 	}
290 	ierr = ISRestoreIndices( proc_is, &is_idx );     CHKERRQ(ierr);
291 	ierr = ISDestroy( &proc_is );                    CHKERRQ(ierr);
292       }
293       ierr = MatDestroy( &adj );                       CHKERRQ(ierr);
294 
295       ierr = ISCreateGeneral( wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc );
296       CHKERRQ(ierr);
297       if( newproc_idx != 0 ) {
298 	ierr = PetscFree( newproc_idx );  CHKERRQ(ierr);
299       }
300     } /* repartitioning */
301     else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
302 
303       PetscInt rfactor,targetPE;
304       /* find factor */
305       if( new_npe == 1 ) rfactor = npe; /* easy */
306       else {
307 	PetscReal best_fact = 0.;
308 	jj = -1;
309 	for( kk = 1 ; kk <= npe ; kk++ ){
310 	  if( npe%kk==0 ) { /* a candidate */
311 	    PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe;
312 	    if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */
313 	    if( fact > best_fact ) {
314 	      best_fact = fact; jj = kk;
315 	    }
316 	  }
317 	}
318 	if( jj != -1 ) rfactor = jj;
319 	else rfactor = 1; /* does this happen .. a prime */
320       }
321       new_npe = npe/rfactor;
322 
323       if( new_npe==nactive ) {
324 	*a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
325 	ierr = PetscFree( counts );  CHKERRQ(ierr);
326 	if (pc_gamg->verbose>0){
327           PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq);
328         }
329 	PetscFunctionReturn(0);
330       }
331 
332       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq);
333       targetPE = mype/rfactor;
334       ierr = ISCreateStride( wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc ); CHKERRQ(ierr);
335 
336       if( stokes ) {
337         ierr = ISCreateStride( wcomm, ncrs_prim*cr_bs, targetPE, 0, &is_eq_newproc_prim ); CHKERRQ(ierr);
338       }
339     } /* end simple 'is_eq_newproc' */
340 
341     /*
342      Create an index set from the is_eq_newproc index set to indicate the mapping TO
343      */
344     ierr = ISPartitioningToNumbering( is_eq_newproc, &is_eq_num ); CHKERRQ(ierr);
345     if( stokes ) {
346       ierr = ISPartitioningToNumbering( is_eq_newproc_prim, &is_eq_num_prim ); CHKERRQ(ierr);
347     }
348     else is_eq_num_prim = is_eq_num;
349     /*
350       Determine how many equations/vertices are assigned to each processor
351      */
352     ierr = ISPartitioningCount( is_eq_newproc, npe, counts ); CHKERRQ(ierr);
353     ncrs_eq_new = counts[mype];
354     ierr = ISDestroy( &is_eq_newproc );                       CHKERRQ(ierr);
355     if( stokes ) {
356       ierr = ISPartitioningCount( is_eq_newproc_prim, npe, counts ); CHKERRQ(ierr);
357       ierr = ISDestroy( &is_eq_newproc_prim );                       CHKERRQ(ierr);
358       ncrs_prim_new = counts[mype]/cr_bs; /* nodes */
359     }
360     else ncrs_prim_new = ncrs_eq_new/cr_bs; /* eqs */
361 
362     ierr = PetscFree( counts );  CHKERRQ(ierr);
363 #if defined PETSC_GAMG_USE_LOG
364     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);   CHKERRQ(ierr);
365 #endif
366 
367     /* move data (for primal equations only) */
368     /* Create a vector to contain the newly ordered element information */
369     ierr = VecCreate( wcomm, &dest_crd );
370     ierr = VecSetSizes( dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE ); CHKERRQ(ierr);
371     ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */
372     /*
373      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
374      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
375      */
376     ierr = PetscMalloc( (ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr);
377     ierr = ISGetIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
378     for(ii=0,jj=0; ii<ncrs_prim ; ii++) {
379       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
380       for( kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
381     }
382     ierr = ISRestoreIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr);
383     ierr = ISCreateGeneral( wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat );
384     CHKERRQ(ierr);
385     ierr = PetscFree( tidx );  CHKERRQ(ierr);
386     /*
387      Create a vector to contain the original vertex information for each element
388      */
389     ierr = VecCreateSeq( PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd ); CHKERRQ(ierr);
390     for( jj=0; jj<ndata_cols ; jj++ ) {
391       const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows;
392       for( ii=0 ; ii<ncrs_prim ; ii++) {
393 	for( kk=0; kk<ndata_rows ; kk++ ) {
394 	  PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
395           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
396 	  ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES );  CHKERRQ(ierr);
397 	}
398       }
399     }
400     ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr);
401     ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr);
402     /*
403       Scatter the element vertex information (still in the original vertex ordering)
404       to the correct processor
405     */
406     ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat);
407     CHKERRQ(ierr);
408     ierr = ISDestroy( &isscat );       CHKERRQ(ierr);
409     ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410     ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
411     ierr = VecScatterDestroy( &vecscat );       CHKERRQ(ierr);
412     ierr = VecDestroy( &src_crd );       CHKERRQ(ierr);
413     /*
414       Put the element vertex data into a new allocation of the gdata->ele
415     */
416     ierr = PetscFree( pc_gamg->data );    CHKERRQ(ierr);
417     ierr = PetscMalloc( node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data );    CHKERRQ(ierr);
418     pc_gamg->data_sz = node_data_sz*ncrs_prim_new;
419     strideNew = ncrs_prim_new*ndata_rows;
420     ierr = VecGetArray( dest_crd, &array );    CHKERRQ(ierr);
421     for( jj=0; jj<ndata_cols ; jj++ ) {
422       for( ii=0 ; ii<ncrs_prim_new ; ii++) {
423 	for( kk=0; kk<ndata_rows ; kk++ ) {
424 	  PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
425 	  pc_gamg->data[ix] = PetscRealPart(array[jx]);
426 	}
427       }
428     }
429     ierr = VecRestoreArray( dest_crd, &array );    CHKERRQ(ierr);
430     ierr = VecDestroy( &dest_crd );    CHKERRQ(ierr);
431 
432     /* move A and P (columns) with new layout */
433 #if defined PETSC_GAMG_USE_LOG
434     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
435 #endif
436 
437     /*
438       Invert for MatGetSubMatrix
439     */
440     ierr = ISInvertPermutation( is_eq_num, ncrs_eq_new, &new_eq_indices ); CHKERRQ(ierr);
441     ierr = ISSort( new_eq_indices ); CHKERRQ(ierr); /* is this needed? */
442     ierr = ISSetBlockSize( new_eq_indices, cr_bs );   CHKERRQ(ierr);
443     if(is_eq_num != is_eq_num_prim) {
444       ierr = ISDestroy( &is_eq_num_prim ); CHKERRQ(ierr); /* could be same as 'is_eq_num' */
445     }
446     ierr = ISDestroy( &is_eq_num ); CHKERRQ(ierr);
447 #if defined PETSC_GAMG_USE_LOG
448     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
449     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
450 #endif
451     /* 'a_Amat_crs' output */
452     {
453       Mat mat;
454       ierr = MatGetSubMatrix( Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat );
455       CHKERRQ(ierr);
456       *a_Amat_crs = mat;
457 
458       if(!PETSC_TRUE){
459         PetscInt cbs, rbs;
460         ierr = MatGetBlockSizes( Cmat, &rbs, &cbs ); CHKERRQ(ierr);
461         PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
462         ierr = MatGetBlockSizes( mat, &rbs, &cbs ); CHKERRQ(ierr);
463         PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",mype,__FUNCT__,rbs,cbs,cr_bs);
464       }
465     }
466     ierr = MatDestroy( &Cmat ); CHKERRQ(ierr);
467 
468 #if defined PETSC_GAMG_USE_LOG
469     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
470 #endif
471     /* prolongator */
472     {
473       IS findices;
474       PetscInt Istart,Iend;
475       Mat Pnew;
476       ierr = MatGetOwnershipRange( Pold, &Istart, &Iend );    CHKERRQ(ierr);
477 #if defined PETSC_GAMG_USE_LOG
478       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
479 #endif
480       ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices);   CHKERRQ(ierr);
481       ierr = ISSetBlockSize(findices,f_bs);   CHKERRQ(ierr);
482       ierr = MatGetSubMatrix( Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew );
483       CHKERRQ(ierr);
484       ierr = ISDestroy( &findices ); CHKERRQ(ierr);
485 
486       if(!PETSC_TRUE){
487         PetscInt cbs, rbs;
488         ierr = MatGetBlockSizes( Pold, &rbs, &cbs ); CHKERRQ(ierr);
489         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
490         ierr = MatGetBlockSizes( Pnew, &rbs, &cbs ); CHKERRQ(ierr);
491         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",mype,__FUNCT__,rbs,cbs);
492       }
493 #if defined PETSC_GAMG_USE_LOG
494       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
495 #endif
496       ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr);
497 
498       /* output - repartitioned */
499       *a_P_inout = Pnew;
500     }
501     ierr = ISDestroy( &new_eq_indices ); CHKERRQ(ierr);
502 
503     *a_nactive_proc = new_npe; /* output */
504   }
505 
506   /* outout matrix data */
507   if( !PETSC_TRUE ) {
508     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
509     if(llev==0) {
510       sprintf(fname,"Cmat_%d.m",llev++);
511       PetscViewerASCIIOpen(wcomm,fname,&viewer);
512       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
513       ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr);
514       ierr = PetscViewerDestroy( &viewer );
515     }
516     sprintf(fname,"Cmat_%d.m",llev++);
517     PetscViewerASCIIOpen(wcomm,fname,&viewer);
518     ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
519     ierr = MatView(Cmat, viewer ); CHKERRQ(ierr);
520     ierr = PetscViewerDestroy( &viewer );
521   }
522 
523   PetscFunctionReturn(0);
524 }
525 
526 /* -------------------------------------------------------------------------- */
527 /*
528    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
529                     by setting data structures and options.
530 
531    Input Parameter:
532 .  pc - the preconditioner context
533 
534    Application Interface Routine: PCSetUp()
535 
536    Notes:
537    The interface routine PCSetUp() is not usually called directly by
538    the user, but instead is called by PCApply() if necessary.
539 */
540 #undef __FUNCT__
541 #define __FUNCT__ "PCSetUp_GAMG"
542 PetscErrorCode PCSetUp_GAMG( PC pc )
543 {
544   PetscErrorCode  ierr;
545   PC_MG           *mg = (PC_MG*)pc->data;
546   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
547   Mat              Pmat = pc->pmat;
548   PetscInt         fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
549   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
550   PetscMPIInt      mype,npe,nactivepe;
551   Mat              Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
552   PetscReal        emaxs[GAMG_MAXLEVELS];
553   IS              *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS];
554   PetscInt         level_bs[GAMG_MAXLEVELS];
555   GAMGKKTMat       kktMatsArr[GAMG_MAXLEVELS];
556   PetscLogDouble   nnz0=0.,nnztot=0.;
557   MatInfo          info;
558   PetscBool        stokes = PETSC_FALSE;
559 
560   PetscFunctionBegin;
561   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
562   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
563   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);
564   if( pc_gamg->setup_count++ > 0 ) {
565     PC_MG_Levels **mglevels = mg->levels;
566     /* just do Galerkin grids */
567     Mat B,dA,dB;
568     assert(pc->setupcalled);
569 
570     if( pc_gamg->Nlevels > 1 ) {
571       /* currently only handle case where mat and pmat are the same on coarser levels */
572       ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr);
573       /* (re)set to get dirty flag */
574       ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
575 
576       for (level=pc_gamg->Nlevels-2; level>-1; level--) {
577         /* the first time through the matrix structure has changed from repartitioning */
578         if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) {
579           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
580           ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
581           mglevels[level]->A = B;
582         }
583         else {
584           ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
585           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
586         }
587         ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
588         dB = B;
589       }
590     }
591 
592     ierr = PCSetUp_MG( pc );CHKERRQ( ierr );
593 
594     /* PCSetUp_MG seems to insists on setting this to GMRES */
595     ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr);
596 
597     PetscFunctionReturn(0);
598   }
599   assert(pc->setupcalled == 0);
600 
601   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
602 
603   ierr = GAMGKKTMatCreate( Pmat, stokes, &kktMatsArr[0] ); CHKERRQ(ierr);
604 
605   if( pc_gamg->data == 0 ) {
606     if( !pc_gamg->createdefaultdata ){
607       SETERRQ(wcomm,PETSC_ERR_LIB,"'createdefaultdata' not set(?) need to support NULL data");
608     }
609     if( stokes ) {
610       SETERRQ(wcomm,PETSC_ERR_LIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
611     }
612     ierr = pc_gamg->createdefaultdata( pc, kktMatsArr[0].A11 ); CHKERRQ(ierr);
613   }
614 
615   /* get basic dims */
616   if( stokes ) {
617     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
618   }
619   else {
620     ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr);
621   }
622 
623   ierr = MatGetSize( Pmat, &M, &qq );CHKERRQ(ierr);
624   if (pc_gamg->verbose) {
625     if(pc_gamg->verbose==1) ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);
626     else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);
627     CHKERRQ(ierr);
628     nnz0 = info.nz_used;
629     nnztot = info.nz_used;
630     PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
631                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
632                 (int)(nnz0/(PetscReal)M),npe);
633   }
634 
635   /* Get A_i and R_i */
636   for ( level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
637         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
638         level++ ){
639     level1 = level + 1;
640 #if defined PETSC_GAMG_USE_LOG
641     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr);
642 #if (defined GAMG_STAGES)
643     ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr );
644 #endif
645 #endif
646     /* deal with Stokes, get sub matrices */
647     if( level > 0 ) {
648       ierr = GAMGKKTMatCreate( Aarr[level], stokes, &kktMatsArr[level] ); CHKERRQ(ierr);
649     }
650     { /* construct prolongator */
651       Mat Gmat;
652       PetscCoarsenData *agg_lists;
653       Mat Prol11,Prol22;
654 
655       level_bs[level] = bs;
656       ierr = pc_gamg->graph( pc,kktMatsArr[level].A11, &Gmat ); CHKERRQ(ierr);
657       ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr);
658       ierr = pc_gamg->prolongator( pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11 ); CHKERRQ(ierr);
659 
660       /* could have failed to create new level */
661       if( Prol11 ){
662         /* get new block size of coarse matrices */
663         ierr = MatGetBlockSizes( Prol11, PETSC_NULL, &bs ); CHKERRQ(ierr);
664 
665         if( stokes ) {
666           if(!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
667           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
668           ierr = pc_gamg->formkktprol( pc, Prol11, kktMatsArr[level].A21, &Prol22 ); CHKERRQ(ierr);
669         }
670 
671         if( pc_gamg->optprol ){
672           /* smooth */
673           ierr = pc_gamg->optprol( pc, kktMatsArr[level].A11, &Prol11 ); CHKERRQ(ierr);
674         }
675 
676         /* remove rows of singleton rows (BCs) */
677 
678         if( stokes ) {
679           IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is};
680           Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 };
681           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1] ); CHKERRQ(ierr);
682         }
683         else {
684           Parr[level1] = Prol11;
685         }
686       }
687       else Parr[level1] = PETSC_NULL;
688 
689       if ( pc_gamg->use_aggs_in_gasm ) {
690         ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level] );  CHKERRQ(ierr);
691       }
692 
693       ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] );  CHKERRQ(ierr);
694 
695       ierr = MatDestroy( &Gmat );      CHKERRQ(ierr);
696       ierr = PetscCDDestroy( agg_lists );  CHKERRQ(ierr);
697     } /* construct prolongator scope */
698 #if defined PETSC_GAMG_USE_LOG
699     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
700 #endif
701     /* cache eigen estimate */
702     if( pc_gamg->emax_id != -1 ){
703       PetscBool flag;
704       ierr = PetscObjectComposedDataGetReal( (PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag );
705       CHKERRQ( ierr );
706       if( !flag ) emaxs[level] = -1.;
707     }
708     else emaxs[level] = -1.;
709     if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
710     if( !Parr[level1] ) {
711       if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);
712       break;
713     }
714 #if defined PETSC_GAMG_USE_LOG
715     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
716 #endif
717 
718     ierr = createLevel( pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
719                         stokes, &Parr[level1], &Aarr[level1], &nactivepe );
720     CHKERRQ(ierr);
721 
722 #if defined PETSC_GAMG_USE_LOG
723     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
724 #endif
725     ierr = MatGetSize( Aarr[level1], &M, &qq );CHKERRQ(ierr);
726 
727     if (pc_gamg->verbose > 0){
728       PetscInt NN = M;
729       if(pc_gamg->verbose==1) {
730         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); CHKERRQ(ierr);
731         ierr = MatGetLocalSize( Aarr[level1], &NN, &qq );
732       }
733       else ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info );
734 
735       CHKERRQ(ierr);
736       nnztot += info.nz_used;
737       PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
738                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
739                   (int)(info.nz_used/(PetscReal)NN), nactivepe );
740       CHKERRQ(ierr);
741     }
742 
743     /* stop if one node -- could pull back for singular problems */
744     if( M/pc_gamg->data_cell_cols < 2 ) {
745       level++;
746       break;
747     }
748 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
749     ierr = PetscLogStagePop(); CHKERRQ( ierr );
750 #endif
751   } /* levels */
752 
753   if( pc_gamg->data ) {
754     ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
755     pc_gamg->data = PETSC_NULL;
756   }
757 
758   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
759   pc_gamg->Nlevels = level + 1;
760   fine_level = level;
761   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
762 
763   /* simple setup */
764   if( !PETSC_TRUE ){
765     PC_MG_Levels **mglevels = mg->levels;
766     for (lidx=0,level=pc_gamg->Nlevels-1;
767          lidx<fine_level;
768          lidx++, level--){
769       ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr);
770       ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr);
771       ierr = MatDestroy( &Parr[level] );  CHKERRQ(ierr);
772       ierr = MatDestroy( &Aarr[level] );  CHKERRQ(ierr);
773     }
774     ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
775 
776     ierr = PCSetUp_MG( pc );  CHKERRQ( ierr );
777   }
778   else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */
779     /* set default smoothers & set operators */
780     for ( lidx = 1, level = pc_gamg->Nlevels-2;
781           lidx <= fine_level;
782           lidx++, level--) {
783       KSP smoother;
784       PC subpc;
785 
786       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
787       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
788 
789       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
790       /* set ops */
791       ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );   CHKERRQ(ierr);
792       ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr);
793 
794       /* create field split PC, get subsmoother */
795       if( stokes ) {
796         KSP *ksps;
797         PetscInt nn;
798         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);   CHKERRQ(ierr);
799         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is); CHKERRQ(ierr);
800         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
801         smoother = ksps[0];
802         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
803         ierr = PetscFree( ksps ); CHKERRQ(ierr);
804       }
805       ierr = GAMGKKTMatDestroy( &kktMatsArr[level] ); CHKERRQ(ierr);
806 
807       /* set defaults */
808       ierr = KSPSetType( smoother, KSPCHEBYSHEV );CHKERRQ(ierr);
809 
810       /* override defaults and command line args (!) */
811       if ( pc_gamg->use_aggs_in_gasm ) {
812         PetscInt sz;
813         IS *is;
814 
815         sz = nASMBlocksArr[level];
816         is = ASMLocalIDsArr[level];
817         ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr);
818         if(sz==0){
819           IS is;
820           PetscInt my0,kk;
821           ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr);
822           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr);
823           ierr = PCGASMSetSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr);
824           ierr = ISDestroy( &is ); CHKERRQ(ierr);
825         }
826         else {
827           PetscInt kk;
828           ierr = PCGASMSetSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr);
829           for(kk=0;kk<sz;kk++){
830             ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr);
831           }
832           ierr = PetscFree( is ); CHKERRQ(ierr);
833         }
834         ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr);
835 
836         ASMLocalIDsArr[level] = PETSC_NULL;
837         nASMBlocksArr[level] = 0;
838         ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr);
839       }
840       else {
841         ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr);
842       }
843     }
844     {
845       /* coarse grid */
846       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
847       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
848       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
849       ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr);
850       ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr);
851       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
852       ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr);
853       ierr = PCSetUp( subpc ); CHKERRQ(ierr);
854       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
855       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
856       ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr);
857     }
858 
859     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
860     ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr);
861     ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr);
862     ierr = PetscOptionsEnd();  CHKERRQ(ierr);
863     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
864 
865     /* create cheby smoothers */
866     for ( lidx = 1, level = pc_gamg->Nlevels-2;
867           lidx <= fine_level;
868           lidx++, level--) {
869       KSP smoother;
870       PetscBool flag;
871       PC subpc;
872 
873       ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr);
874       ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
875 
876       /* create field split PC, get subsmoother */
877       if( stokes ) {
878         KSP *ksps;
879         PetscInt nn;
880         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr);
881         smoother = ksps[0];
882         ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr);
883         ierr = PetscFree( ksps );  CHKERRQ(ierr);
884       }
885 
886       /* do my own cheby */
887       ierr = PetscObjectTypeCompare( (PetscObject)smoother, KSPCHEBYSHEV, &flag ); CHKERRQ(ierr);
888       if( flag ) {
889         PetscReal emax, emin;
890         ierr = PetscObjectTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr);
891         if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
892         else{ /* eigen estimate 'emax' */
893           KSP eksp; Mat Lmat = Aarr[level];
894           Vec bb, xx;
895 
896           ierr = MatGetVecs( Lmat, &bb, 0 );         CHKERRQ(ierr);
897           ierr = MatGetVecs( Lmat, &xx, 0 );         CHKERRQ(ierr);
898           {
899             PetscRandom    rctx;
900             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
901             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
902             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
903             ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr);
904           }
905 
906           if( removedEqs[level] ) {
907             /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */
908             PetscScalar *zeros;
909             PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level];
910             const PetscInt *idx;
911             ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr);
912             ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr);
913             for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.;
914             ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr);
915             ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr);
916             for(ii=0;ii<sz;ii++) {
917               for(jj=0;jj<bs;jj++) {
918                 idx_bs[ii] = bs*idx[ii]+jj;
919               }
920             }
921             ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr);
922             if( sz > 0 ) {
923               ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES );  CHKERRQ(ierr);
924             }
925             ierr = PetscFree( idx_bs );  CHKERRQ(ierr);
926             ierr = PetscFree( zeros );  CHKERRQ(ierr);
927             ierr = VecAssemblyBegin(bb); CHKERRQ(ierr);
928             ierr = VecAssemblyEnd(bb); CHKERRQ(ierr);
929           }
930           ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr);
931           ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 );
932           CHKERRQ(ierr);
933           ierr = KSPSetNormType( eksp, KSP_NORM_NONE );                 CHKERRQ(ierr);
934           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
935           ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_");         CHKERRQ(ierr);
936           ierr = KSPSetFromOptions( eksp );    CHKERRQ(ierr);
937 
938           ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr);
939           ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr );
940           ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr);
941 
942           /* set PC type to be same as smoother */
943           ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr );
944 
945           /* solve - keep stuff out of logging */
946           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
947           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
948           ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr);
949           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
950           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
951 
952           ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr);
953 
954           ierr = VecDestroy( &xx );       CHKERRQ(ierr);
955           ierr = VecDestroy( &bb );       CHKERRQ(ierr);
956           ierr = KSPDestroy( &eksp );       CHKERRQ(ierr);
957 
958           if( pc_gamg->verbose > 0 ) {
959             PetscInt N1, tt;
960             ierr = MatGetSize( Aarr[level], &N1, &tt );         CHKERRQ(ierr);
961             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
962           }
963         }
964         {
965           PetscInt N1, N0;
966           ierr = MatGetSize( Aarr[level], &N1, PETSC_NULL );         CHKERRQ(ierr);
967           ierr = MatGetSize( Aarr[level+1], &N0, PETSC_NULL );       CHKERRQ(ierr);
968           /* heuristic - is this crap? */
969           emin = 1.*emax/((PetscReal)N1/(PetscReal)N0);
970           emax *= 1.05;
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   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1357 
1358   PetscFunctionBegin;
1359   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1360   {
1361     /* -pc_gamg_verbose */
1362     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
1363                            "none", pc_gamg->verbose,
1364                            &pc_gamg->verbose, PETSC_NULL );
1365     CHKERRQ(ierr);
1366 
1367     /* -pc_gamg_repartition */
1368     ierr = PetscOptionsBool("-pc_gamg_repartition",
1369                             "Repartion coarse grids (false)",
1370                             "PCGAMGRepartitioning",
1371                             pc_gamg->repart,
1372                             &pc_gamg->repart,
1373                             &flag);
1374     CHKERRQ(ierr);
1375 
1376     /* -pc_gamg_use_agg_gasm */
1377     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1378                             "Use aggregation agragates for GASM smoother (false)",
1379                             "PCGAMGUseASMAggs",
1380                             pc_gamg->use_aggs_in_gasm,
1381                             &pc_gamg->use_aggs_in_gasm,
1382                             &flag);
1383     CHKERRQ(ierr);
1384 
1385     /* -pc_gamg_process_eq_limit */
1386     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1387                            "Limit (goal) on number of equations per process on coarse grids",
1388                            "PCGAMGSetProcEqLim",
1389                            pc_gamg->min_eq_proc,
1390                            &pc_gamg->min_eq_proc,
1391                            &flag );
1392     CHKERRQ(ierr);
1393 
1394     /* -pc_gamg_coarse_eq_limit */
1395     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1396                            "Limit on number of equations for the coarse grid",
1397                            "PCGAMGSetCoarseEqLim",
1398                            pc_gamg->coarse_eq_limit,
1399                            &pc_gamg->coarse_eq_limit,
1400                            &flag );
1401     CHKERRQ(ierr);
1402 
1403     /* -pc_gamg_threshold */
1404     ierr = PetscOptionsReal("-pc_gamg_threshold",
1405                             "Relative threshold to use for dropping edges in aggregation graph",
1406                             "PCGAMGSetThreshold",
1407                             pc_gamg->threshold,
1408                             &pc_gamg->threshold,
1409                             &flag );
1410     CHKERRQ(ierr);
1411     if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);
1412 
1413     ierr = PetscOptionsInt("-pc_mg_levels",
1414                            "Set number of MG levels",
1415                            "PCGAMGSetNlevels",
1416                            pc_gamg->Nlevels,
1417                            &pc_gamg->Nlevels,
1418                            &flag );
1419   }
1420   ierr = PetscOptionsTail();CHKERRQ(ierr);
1421 
1422   PetscFunctionReturn(0);
1423 }
1424 
1425 /* -------------------------------------------------------------------------- */
1426 /*MC
1427      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
1428        - This is the entry point to GAMG, registered in pcregis.c
1429 
1430    Options Database Keys:
1431    Multigrid options(inherited)
1432 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1433 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1434 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1435 -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
1436 
1437   Level: intermediate
1438 
1439   Concepts: multigrid
1440 
1441 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1442            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1443            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1444            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1445            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1446 M*/
1447 EXTERN_C_BEGIN
1448 #undef __FUNCT__
1449 #define __FUNCT__ "PCCreate_GAMG"
1450 PetscErrorCode  PCCreate_GAMG( PC pc )
1451 {
1452   PetscErrorCode  ierr;
1453   PC_GAMG         *pc_gamg;
1454   PC_MG           *mg;
1455 #if defined PETSC_GAMG_USE_LOG
1456   static long count = 0;
1457 #endif
1458 
1459   PetscFunctionBegin;
1460 
1461   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1462   ierr = PCSetType( pc, PCMG );  CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
1463   ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr);
1464 
1465   /* create a supporting struct and attach it to pc */
1466   ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr);
1467   mg = (PC_MG*)pc->data;
1468   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
1469   mg->innerctx = pc_gamg;
1470 
1471   pc_gamg->setup_count = 0;
1472   /* these should be in subctx but repartitioning needs simple arrays */
1473   pc_gamg->data_sz = 0;
1474   pc_gamg->data = 0;
1475 
1476   /* register AMG type */
1477   if( !GAMGList ){
1478     ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1479     ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
1480   }
1481 
1482   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1483   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1484   pc->ops->setup          = PCSetUp_GAMG;
1485   pc->ops->reset          = PCReset_GAMG;
1486   pc->ops->destroy        = PCDestroy_GAMG;
1487 
1488   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1489 					    "PCGAMGSetProcEqLim_C",
1490 					    "PCGAMGSetProcEqLim_GAMG",
1491 					    PCGAMGSetProcEqLim_GAMG);
1492   CHKERRQ(ierr);
1493 
1494   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1495 					    "PCGAMGSetCoarseEqLim_C",
1496 					    "PCGAMGSetCoarseEqLim_GAMG",
1497 					    PCGAMGSetCoarseEqLim_GAMG);
1498   CHKERRQ(ierr);
1499 
1500   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1501 					    "PCGAMGSetRepartitioning_C",
1502 					    "PCGAMGSetRepartitioning_GAMG",
1503 					    PCGAMGSetRepartitioning_GAMG);
1504   CHKERRQ(ierr);
1505 
1506   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1507 					    "PCGAMGSetUseASMAggs_C",
1508 					    "PCGAMGSetUseASMAggs_GAMG",
1509 					    PCGAMGSetUseASMAggs_GAMG);
1510   CHKERRQ(ierr);
1511 
1512   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1513 					    "PCGAMGSetThreshold_C",
1514 					    "PCGAMGSetThreshold_GAMG",
1515 					    PCGAMGSetThreshold_GAMG);
1516   CHKERRQ(ierr);
1517 
1518   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
1519 					    "PCGAMGSetType_C",
1520 					    "PCGAMGSetType_GAMG",
1521 					    PCGAMGSetType_GAMG);
1522   CHKERRQ(ierr);
1523 
1524   pc_gamg->repart = PETSC_FALSE;
1525   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1526   pc_gamg->min_eq_proc = 100;
1527   pc_gamg->coarse_eq_limit = 800;
1528   pc_gamg->threshold = 0.001;
1529   pc_gamg->Nlevels = GAMG_MAXLEVELS;
1530   pc_gamg->verbose = 0;
1531   pc_gamg->emax_id = -1;
1532 
1533   /* private events */
1534 #if defined PETSC_GAMG_USE_LOG
1535   if( count++ == 0 ) {
1536     PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1537     PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1538     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1539     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1540     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1541     PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1542     PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1543     PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1544     PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1545     PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1546     PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1547     PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1548     PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1549     PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1550     PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1551     PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1552     PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1553 
1554     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1555     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1556     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1557     /* create timer stages */
1558 #if defined GAMG_STAGES
1559     {
1560       char str[32];
1561       sprintf(str,"MG Level %d (finest)",0);
1562       PetscLogStageRegister(str, &gamg_stages[0]);
1563       PetscInt lidx;
1564       for (lidx=1;lidx<9;lidx++){
1565 	sprintf(str,"MG Level %d",lidx);
1566 	PetscLogStageRegister(str, &gamg_stages[lidx]);
1567       }
1568     }
1569 #endif
1570   }
1571 #endif
1572   /* general events */
1573 #if defined PETSC_USE_LOG
1574   PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);
1575   PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);
1576   PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1577   PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1578   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1579   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1580   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);
1581   PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);
1582 #endif
1583 
1584   /* instantiate derived type */
1585   ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr);
1586   {
1587     char tname[256] = GAMGAGG;
1588     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",
1589                             GAMGList, tname, tname, sizeof(tname), PETSC_NULL );
1590     CHKERRQ(ierr);
1591     ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr);
1592   }
1593   ierr = PetscOptionsTail();   CHKERRQ(ierr);
1594 
1595   PetscFunctionReturn(0);
1596 }
1597 EXTERN_C_END
1598