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