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