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