xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 5c8f6a953e7ed1c81f507d64aebddb11080b60e9)
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) SETERRQ(wcomm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
611       if (stokes) SETERRQ(wcomm,PETSC_ERR_PLIB,"Need data (eg, PCSetCoordinates) for Stokes problems");
612       ierr = pc_gamg->createdefaultdata(pc, kktMatsArr[0].A11);CHKERRQ(ierr);
613     }
614   }
615 
616   /* cache original data for reuse */
617   if (!pc_gamg->orig_data && redo_mesh_setup) {
618     ierr = PetscMalloc(pc_gamg->data_sz*sizeof(PetscReal), &pc_gamg->orig_data);CHKERRQ(ierr);
619     for (qq=0;qq<pc_gamg->data_sz;qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
620     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
621     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
622   }
623 
624   /* get basic dims */
625   if (stokes) {
626     bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */
627   } else {
628     ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
629   }
630 
631   ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr);
632   if (pc_gamg->verbose) {
633     PetscInt NN = M;
634     if (pc_gamg->verbose==1) {
635       ierr =  MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr);
636       ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr);
637     } else {
638       ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
639     }
640     nnz0 = info.nz_used;
641     nnztot = info.nz_used;
642     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",
643                 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
644                 (int)(nnz0/(PetscReal)NN),npe);CHKERRQ(ierr);
645   }
646 
647   /* Get A_i and R_i */
648   for (level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */
649         level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */
650         level++) {
651     level1 = level + 1;
652 #if defined PETSC_GAMG_USE_LOG
653     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
654 #if (defined GAMG_STAGES)
655     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
656 #endif
657 #endif
658     /* deal with Stokes, get sub matrices */
659     if (level > 0) {
660       ierr = GAMGKKTMatCreate(Aarr[level], stokes, &kktMatsArr[level]);CHKERRQ(ierr);
661     }
662     { /* construct prolongator */
663       Mat Gmat;
664       PetscCoarsenData *agg_lists;
665       Mat Prol11,Prol22;
666 
667       ierr = pc_gamg->graph(pc,kktMatsArr[level].A11, &Gmat);CHKERRQ(ierr);
668       ierr = pc_gamg->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
669       ierr = pc_gamg->prolongator(pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11);CHKERRQ(ierr);
670 
671       /* could have failed to create new level */
672       if (Prol11) {
673         /* get new block size of coarse matrices */
674         ierr = MatGetBlockSizes(Prol11, PETSC_NULL, &bs);CHKERRQ(ierr);
675 
676         if (stokes) {
677           if (!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method.");
678           /* R A12 == (T = A21 P)';  G = T' T; coarsen G; form plain agg with G */
679           ierr = pc_gamg->formkktprol(pc, Prol11, kktMatsArr[level].A21, &Prol22);CHKERRQ(ierr);
680         }
681 
682         if (pc_gamg->optprol) {
683           /* smooth */
684           ierr = pc_gamg->optprol(pc, kktMatsArr[level].A11, &Prol11);CHKERRQ(ierr);
685         }
686 
687         if (stokes) {
688           IS is_row[2];
689           Mat a[4];
690 
691           is_row[0] = kktMatsArr[level].prim_is; is_row[1] = kktMatsArr[level].constr_is;
692           a[0] = Prol11; a[1] = PETSC_NULL; a[2] = PETSC_NULL; a[3] = Prol22;
693           ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1]);CHKERRQ(ierr);
694         } else {
695           Parr[level1] = Prol11;
696         }
697       } else Parr[level1] = PETSC_NULL;
698 
699       if (pc_gamg->use_aggs_in_gasm) {
700         ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
701       }
702 
703       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
704       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
705     } /* construct prolongator scope */
706 #if defined PETSC_GAMG_USE_LOG
707     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
708 #endif
709     /* cache eigen estimate */
710     if (pc_gamg->emax_id != -1) {
711       PetscBool flag;
712       ierr = PetscObjectComposedDataGetReal((PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr);
713       if (!flag) emaxs[level] = -1.;
714     } else emaxs[level] = -1.;
715     if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
716     if (!Parr[level1]) {
717       if (pc_gamg->verbose) {
718         ierr =  PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level);CHKERRQ(ierr);
719       }
720       break;
721     }
722 #if defined PETSC_GAMG_USE_LOG
723     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
724 #endif
725 
726     ierr = createLevel(pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2),
727                         stokes, &Parr[level1], &Aarr[level1], &nactivepe);CHKERRQ(ierr);
728 
729 #if defined PETSC_GAMG_USE_LOG
730     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
731 #endif
732     ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr);
733 
734     if (pc_gamg->verbose > 0) {
735       PetscInt NN = M;
736       if (pc_gamg->verbose==1) {
737         ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr);
738         ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr);
739       } else {
740         ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
741       }
742 
743       nnztot += info.nz_used;
744       ierr = PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",
745                   mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols,
746                   (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr);
747     }
748 
749     /* stop if one node -- could pull back for singular problems */
750     if (M/pc_gamg->data_cell_cols < 2) {
751       level++;
752       break;
753     }
754 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
755     ierr = PetscLogStagePop();CHKERRQ(ierr);
756 #endif
757   } /* levels */
758 
759   if (pc_gamg->data) {
760     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
761     pc_gamg->data = PETSC_NULL;
762   }
763 
764   if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0);
765   pc_gamg->Nlevels = level + 1;
766   fine_level = level;
767   ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr);
768 
769   /* simple setup */
770   if (!PETSC_TRUE) {
771     PC_MG_Levels **mglevels = mg->levels;
772     for (lidx=0,level=pc_gamg->Nlevels-1;
773          lidx<fine_level;
774          lidx++, level--) {
775       ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr);
776       ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
777       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
778       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
779     }
780     ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
781 
782     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
783   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
784     /* set default smoothers & set operators */
785     for (lidx = 1, level = pc_gamg->Nlevels-2;
786           lidx <= fine_level;
787           lidx++, level--) {
788       KSP smoother;
789       PC subpc;
790 
791       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
792       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
793 
794       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
795       /* set ops */
796       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
797       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
798 
799       /* create field split PC, get subsmoother */
800       if (stokes) {
801         KSP *ksps;
802         PetscInt nn;
803         ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is);CHKERRQ(ierr);
804         ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is);CHKERRQ(ierr);
805         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
806         smoother = ksps[0];
807         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
808         ierr = PetscFree(ksps);CHKERRQ(ierr);
809       }
810       ierr = GAMGKKTMatDestroy(&kktMatsArr[level]);CHKERRQ(ierr);
811 
812       /* set defaults */
813       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
814 
815       /* override defaults and command line args (!) */
816       if (pc_gamg->use_aggs_in_gasm) {
817         PetscInt sz;
818         IS *is;
819 
820         sz = nASMBlocksArr[level];
821         is = ASMLocalIDsArr[level];
822         ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr);
823         if (sz==0) {
824           IS is;
825           PetscInt my0,kk;
826           ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr);
827           ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
828           ierr = PCGASMSetSubdomains(subpc, 1, &is, PETSC_NULL);CHKERRQ(ierr);
829           ierr = ISDestroy(&is);CHKERRQ(ierr);
830         } else {
831           PetscInt kk;
832           ierr = PCGASMSetSubdomains(subpc, sz, is, PETSC_NULL);CHKERRQ(ierr);
833           for (kk=0;kk<sz;kk++) {
834             ierr = ISDestroy(&is[kk]);CHKERRQ(ierr);
835           }
836           ierr = PetscFree(is);CHKERRQ(ierr);
837         }
838         ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr);
839 
840         ASMLocalIDsArr[level] = PETSC_NULL;
841         nASMBlocksArr[level] = 0;
842         ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr);
843       } else {
844         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
845       }
846     }
847     {
848       /* coarse grid */
849       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
850       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
851       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
852       ierr = KSPSetOperators(smoother, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
853       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
854       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
855       ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
856       ierr = PCSetUp(subpc);CHKERRQ(ierr);
857       ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);      assert(ii==1);
858       ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
859       ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
860       ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
861     }
862 
863     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
864     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
865     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr);
866     ierr = PetscOptionsEnd();CHKERRQ(ierr);
867     if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used.");
868 
869     /* create cheby smoothers */
870     for (lidx = 1, level = pc_gamg->Nlevels-2;
871           lidx <= fine_level;
872           lidx++, level--) {
873       KSP smoother;
874       PetscBool flag;
875       PC subpc;
876 
877       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
878       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
879 
880       /* create field split PC, get subsmoother */
881       if (stokes) {
882         KSP *ksps;
883         PetscInt nn;
884         ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr);
885         smoother = ksps[0];
886         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
887         ierr = PetscFree(ksps);CHKERRQ(ierr);
888       }
889 
890       /* do my own cheby */
891       ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr);
892       if (flag) {
893         PetscReal emax, emin;
894         ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr);
895         if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */
896         else { /* eigen estimate 'emax' */
897           KSP eksp;
898           Mat Lmat = Aarr[level];
899           Vec bb, xx;
900 
901           ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr);
902           ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr);
903           {
904             PetscRandom    rctx;
905             ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr);
906             ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
907             ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr);
908             ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
909           }
910 
911           /* zeroing out BC rows -- needed for crazy matrices */
912           {
913             PetscInt Istart,Iend,ncols,jj,Ii;
914             PetscScalar zero = 0.0;
915             ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr);
916             for (Ii = Istart, jj = 0 ; Ii < Iend ; Ii++, jj++) {
917               ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
918               if (ncols <= 1) {
919                 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr);
920               }
921               ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr);
922             }
923             ierr = VecAssemblyBegin(bb);CHKERRQ(ierr);
924             ierr = VecAssemblyEnd(bb);CHKERRQ(ierr);
925           }
926 
927           ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr);
928           ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr);
929           ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
930           ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
931           ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
932           ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
933 
934           ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
935           ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
936           ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
937 
938           /* set PC type to be same as smoother */
939           ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr);
940 
941           /* solve - keep stuff out of logging */
942           ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
943           ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
944           ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
945           ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
946           ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
947 
948           ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
949 
950           ierr = VecDestroy(&xx);CHKERRQ(ierr);
951           ierr = VecDestroy(&bb);CHKERRQ(ierr);
952           ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
953 
954           if (pc_gamg->verbose > 0) {
955             PetscInt N1, tt;
956             ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr);
957             PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1);
958           }
959         }
960         {
961           PetscInt N1, N0;
962           ierr = MatGetSize(Aarr[level], &N1, PETSC_NULL);CHKERRQ(ierr);
963           ierr = MatGetSize(Aarr[level+1], &N0, PETSC_NULL);CHKERRQ(ierr);
964           /* heuristic - is this crap? */
965           /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */
966           emin = emax * pc_gamg->eigtarget[0];
967           emax *= pc_gamg->eigtarget[1];
968         }
969         ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr);
970       } /* setup checby flag */
971     } /* non-coarse levels */
972 
973     /* clean up */
974     for (level=1;level<pc_gamg->Nlevels;level++) {
975       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
976       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
977     }
978 
979     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
980 
981     if (PETSC_TRUE) {
982       KSP smoother;  /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */
983       ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
984       ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
985     }
986   } else {
987     KSP smoother;
988     if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__);
989     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
990     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr);
991     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
992     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
993   }
994 
995   PetscFunctionReturn(0);
996 }
997 
998 /* ------------------------------------------------------------------------- */
999 /*
1000  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
1001    that was created with PCCreate_GAMG().
1002 
1003    Input Parameter:
1004 .  pc - the preconditioner context
1005 
1006    Application Interface Routine: PCDestroy()
1007 */
1008 #undef __FUNCT__
1009 #define __FUNCT__ "PCDestroy_GAMG"
1010 PetscErrorCode PCDestroy_GAMG(PC pc)
1011 {
1012   PetscErrorCode  ierr;
1013   PC_MG           *mg = (PC_MG*)pc->data;
1014   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
1015 
1016   PetscFunctionBegin;
1017   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
1018   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
1019   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 
1024 #undef __FUNCT__
1025 #define __FUNCT__ "PCGAMGSetProcEqLim"
1026 /*@
1027    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via
1028    processor reduction.
1029 
1030    Not Collective on PC
1031 
1032    Input Parameters:
1033 .  pc - the preconditioner context
1034 
1035 
1036    Options Database Key:
1037 .  -pc_gamg_process_eq_limit
1038 
1039    Level: intermediate
1040 
1041    Concepts: Unstructured multrigrid preconditioner
1042 
1043 .seealso: ()
1044 @*/
1045 PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
1046 {
1047   PetscErrorCode ierr;
1048 
1049   PetscFunctionBegin;
1050   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1051   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 EXTERN_C_BEGIN
1056 #undef __FUNCT__
1057 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG"
1058 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1059 {
1060   PC_MG           *mg = (PC_MG*)pc->data;
1061   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1062 
1063   PetscFunctionBegin;
1064   if (n>0) pc_gamg->min_eq_proc = n;
1065   PetscFunctionReturn(0);
1066 }
1067 EXTERN_C_END
1068 
1069 #undef __FUNCT__
1070 #define __FUNCT__ "PCGAMGSetCoarseEqLim"
1071 /*@
1072    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.
1073 
1074  Collective on PC
1075 
1076    Input Parameters:
1077 .  pc - the preconditioner context
1078 
1079 
1080    Options Database Key:
1081 .  -pc_gamg_coarse_eq_limit
1082 
1083    Level: intermediate
1084 
1085    Concepts: Unstructured multrigrid preconditioner
1086 
1087 .seealso: ()
1088  @*/
1089 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1090 {
1091   PetscErrorCode ierr;
1092 
1093   PetscFunctionBegin;
1094   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1095   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 EXTERN_C_BEGIN
1100 #undef __FUNCT__
1101 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG"
1102 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1103 {
1104   PC_MG           *mg = (PC_MG*)pc->data;
1105   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1106 
1107   PetscFunctionBegin;
1108   if (n>0) pc_gamg->coarse_eq_limit = n;
1109   PetscFunctionReturn(0);
1110 }
1111 EXTERN_C_END
1112 
1113 #undef __FUNCT__
1114 #define __FUNCT__ "PCGAMGSetRepartitioning"
1115 /*@
1116    PCGAMGSetRepartitioning - Repartition the coarse grids
1117 
1118    Collective on PC
1119 
1120    Input Parameters:
1121 .  pc - the preconditioner context
1122 
1123 
1124    Options Database Key:
1125 .  -pc_gamg_repartition
1126 
1127    Level: intermediate
1128 
1129    Concepts: Unstructured multrigrid preconditioner
1130 
1131 .seealso: ()
1132 @*/
1133 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
1134 {
1135   PetscErrorCode ierr;
1136 
1137   PetscFunctionBegin;
1138   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1139   ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1140   PetscFunctionReturn(0);
1141 }
1142 
1143 EXTERN_C_BEGIN
1144 #undef __FUNCT__
1145 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG"
1146 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
1147 {
1148   PC_MG           *mg = (PC_MG*)pc->data;
1149   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1150 
1151   PetscFunctionBegin;
1152   pc_gamg->repart = n;
1153   PetscFunctionReturn(0);
1154 }
1155 EXTERN_C_END
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "PCGAMGSetReuseProl"
1159 /*@
1160    PCGAMGSetReuseProl - Reuse prlongation
1161 
1162    Collective on PC
1163 
1164    Input Parameters:
1165 .  pc - the preconditioner context
1166 
1167 
1168    Options Database Key:
1169 .  -pc_gamg_reuse_interpolation
1170 
1171    Level: intermediate
1172 
1173    Concepts: Unstructured multrigrid preconditioner
1174 
1175 .seealso: ()
1176 @*/
1177 PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n)
1178 {
1179   PetscErrorCode ierr;
1180 
1181   PetscFunctionBegin;
1182   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1183   ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 EXTERN_C_BEGIN
1188 #undef __FUNCT__
1189 #define __FUNCT__ "PCGAMGSetReuseProl_GAMG"
1190 PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n)
1191 {
1192   PC_MG           *mg = (PC_MG*)pc->data;
1193   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1194 
1195   PetscFunctionBegin;
1196   pc_gamg->reuse_prol = n;
1197   PetscFunctionReturn(0);
1198 }
1199 EXTERN_C_END
1200 
1201 #undef __FUNCT__
1202 #define __FUNCT__ "PCGAMGSetUseASMAggs"
1203 /*@
1204    PCGAMGSetUseASMAggs -
1205 
1206    Collective on PC
1207 
1208    Input Parameters:
1209 .  pc - the preconditioner context
1210 
1211 
1212    Options Database Key:
1213 .  -pc_gamg_use_agg_gasm
1214 
1215    Level: intermediate
1216 
1217    Concepts: Unstructured multrigrid preconditioner
1218 
1219 .seealso: ()
1220 @*/
1221 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
1222 {
1223   PetscErrorCode ierr;
1224 
1225   PetscFunctionBegin;
1226   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1227   ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1228   PetscFunctionReturn(0);
1229 }
1230 
1231 EXTERN_C_BEGIN
1232 #undef __FUNCT__
1233 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG"
1234 PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
1235 {
1236   PC_MG           *mg = (PC_MG*)pc->data;
1237   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1238 
1239   PetscFunctionBegin;
1240   pc_gamg->use_aggs_in_gasm = n;
1241   PetscFunctionReturn(0);
1242 }
1243 EXTERN_C_END
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "PCGAMGSetNlevels"
1247 /*@
1248    PCGAMGSetNlevels -
1249 
1250    Not collective on PC
1251 
1252    Input Parameters:
1253 .  pc - the preconditioner context
1254 
1255 
1256    Options Database Key:
1257 .  -pc_mg_levels
1258 
1259    Level: intermediate
1260 
1261    Concepts: Unstructured multrigrid preconditioner
1262 
1263 .seealso: ()
1264 @*/
1265 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1266 {
1267   PetscErrorCode ierr;
1268 
1269   PetscFunctionBegin;
1270   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1271   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
1272   PetscFunctionReturn(0);
1273 }
1274 
1275 EXTERN_C_BEGIN
1276 #undef __FUNCT__
1277 #define __FUNCT__ "PCGAMGSetNlevels_GAMG"
1278 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1279 {
1280   PC_MG           *mg = (PC_MG*)pc->data;
1281   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1282 
1283   PetscFunctionBegin;
1284   pc_gamg->Nlevels = n;
1285   PetscFunctionReturn(0);
1286 }
1287 EXTERN_C_END
1288 
1289 #undef __FUNCT__
1290 #define __FUNCT__ "PCGAMGSetThreshold"
1291 /*@
1292    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1293 
1294    Not collective on PC
1295 
1296    Input Parameters:
1297 .  pc - the preconditioner context
1298 
1299 
1300    Options Database Key:
1301 .  -pc_gamg_threshold
1302 
1303    Level: intermediate
1304 
1305    Concepts: Unstructured multrigrid preconditioner
1306 
1307 .seealso: ()
1308 @*/
1309 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
1310 {
1311   PetscErrorCode ierr;
1312 
1313   PetscFunctionBegin;
1314   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1315   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr);
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 EXTERN_C_BEGIN
1320 #undef __FUNCT__
1321 #define __FUNCT__ "PCGAMGSetThreshold_GAMG"
1322 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
1323 {
1324   PC_MG           *mg = (PC_MG*)pc->data;
1325   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1326 
1327   PetscFunctionBegin;
1328   pc_gamg->threshold = n;
1329   PetscFunctionReturn(0);
1330 }
1331 EXTERN_C_END
1332 
1333 #undef __FUNCT__
1334 #define __FUNCT__ "PCGAMGSetType"
1335 /*@
1336    PCGAMGSetType - Set solution method - calls sub create method
1337 
1338    Collective on PC
1339 
1340    Input Parameters:
1341 .  pc - the preconditioner context
1342 
1343 
1344    Options Database Key:
1345 .  -pc_gamg_type
1346 
1347    Level: intermediate
1348 
1349    Concepts: Unstructured multrigrid preconditioner
1350 
1351 .seealso: ()
1352 @*/
1353 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1354 {
1355   PetscErrorCode ierr;
1356 
1357   PetscFunctionBegin;
1358   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1359   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1360   PetscFunctionReturn(0);
1361 }
1362 
1363 EXTERN_C_BEGIN
1364 #undef __FUNCT__
1365 #define __FUNCT__ "PCGAMGSetType_GAMG"
1366 PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1367 {
1368   PetscErrorCode ierr,(*r)(PC);
1369 
1370   PetscFunctionBegin;
1371   ierr = PetscFunctionListFind(((PetscObject)pc)->comm,GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr);
1372   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1373   ierr = (*r)(pc);CHKERRQ(ierr);
1374   PetscFunctionReturn(0);
1375 }
1376 EXTERN_C_END
1377 
1378 #undef __FUNCT__
1379 #define __FUNCT__ "PCSetFromOptions_GAMG"
1380 PetscErrorCode PCSetFromOptions_GAMG(PC pc)
1381 {
1382   PetscErrorCode  ierr;
1383   PC_MG           *mg = (PC_MG*)pc->data;
1384   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1385   PetscBool        flag;
1386   PetscInt         two = 2;
1387   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
1388 
1389   PetscFunctionBegin;
1390   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1391   {
1392     /* -pc_gamg_verbose */
1393     ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG",
1394                            "none", pc_gamg->verbose,
1395                            &pc_gamg->verbose, PETSC_NULL);CHKERRQ(ierr);
1396     /* -pc_gamg_repartition */
1397     ierr = PetscOptionsBool("-pc_gamg_repartition",
1398                             "Repartion coarse grids (false)",
1399                             "PCGAMGRepartitioning",
1400                             pc_gamg->repart,
1401                             &pc_gamg->repart,
1402                             &flag);CHKERRQ(ierr);
1403     /* -pc_gamg_reuse_interpolation */
1404     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation",
1405                             "Reuse prolongation operator (true)",
1406                             "PCGAMGReuseProl",
1407                             pc_gamg->reuse_prol,
1408                             &pc_gamg->reuse_prol,
1409                             &flag);CHKERRQ(ierr);
1410     /* -pc_gamg_use_agg_gasm */
1411     ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm",
1412                             "Use aggregation agragates for GASM smoother (false)",
1413                             "PCGAMGUseASMAggs",
1414                             pc_gamg->use_aggs_in_gasm,
1415                             &pc_gamg->use_aggs_in_gasm,
1416                             &flag);CHKERRQ(ierr);
1417     /* -pc_gamg_process_eq_limit */
1418     ierr = PetscOptionsInt("-pc_gamg_process_eq_limit",
1419                            "Limit (goal) on number of equations per process on coarse grids",
1420                            "PCGAMGSetProcEqLim",
1421                            pc_gamg->min_eq_proc,
1422                            &pc_gamg->min_eq_proc,
1423                            &flag);CHKERRQ(ierr);
1424     /* -pc_gamg_coarse_eq_limit */
1425     ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit",
1426                            "Limit on number of equations for the coarse grid",
1427                            "PCGAMGSetCoarseEqLim",
1428                            pc_gamg->coarse_eq_limit,
1429                            &pc_gamg->coarse_eq_limit,
1430                            &flag);CHKERRQ(ierr);
1431     /* -pc_gamg_threshold */
1432     ierr = PetscOptionsReal("-pc_gamg_threshold",
1433                             "Relative threshold to use for dropping edges in aggregation graph",
1434                             "PCGAMGSetThreshold",
1435                             pc_gamg->threshold,
1436                             &pc_gamg->threshold,
1437                             &flag);CHKERRQ(ierr);
1438     if (flag && pc_gamg->verbose) {
1439       ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr);
1440     }
1441 
1442     ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,PETSC_NULL);CHKERRQ(ierr);
1443     ierr = PetscOptionsInt("-pc_mg_levels",
1444                            "Set number of MG levels",
1445                            "PCGAMGSetNlevels",
1446                            pc_gamg->Nlevels,
1447                            &pc_gamg->Nlevels,
1448                            &flag);CHKERRQ(ierr);
1449   }
1450   ierr = PetscOptionsTail();CHKERRQ(ierr);
1451 
1452   PetscFunctionReturn(0);
1453 }
1454 
1455 /* -------------------------------------------------------------------------- */
1456 /*MC
1457      PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework.
1458        - This is the entry point to GAMG, registered in pcregis.c
1459 
1460    Options Database Keys:
1461    Multigrid options(inherited)
1462 +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType)
1463 .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1464 .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1465 -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1466 
1467   Level: intermediate
1468 
1469   Concepts: multigrid
1470 
1471 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1472            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1473            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1474            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1475            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1476 M*/
1477 EXTERN_C_BEGIN
1478 #undef __FUNCT__
1479 #define __FUNCT__ "PCCreate_GAMG"
1480 PetscErrorCode  PCCreate_GAMG(PC pc)
1481 {
1482   PetscErrorCode  ierr;
1483   PC_GAMG         *pc_gamg;
1484   PC_MG           *mg;
1485 #if defined PETSC_GAMG_USE_LOG
1486   static long count = 0;
1487 #endif
1488 
1489   PetscFunctionBegin;
1490   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1491   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
1492   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
1493 
1494   /* create a supporting struct and attach it to pc */
1495   ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr);
1496   mg = (PC_MG*)pc->data;
1497   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
1498   mg->innerctx = pc_gamg;
1499 
1500   pc_gamg->setup_count = 0;
1501   /* these should be in subctx but repartitioning needs simple arrays */
1502   pc_gamg->data_sz = 0;
1503   pc_gamg->data = 0;
1504 
1505   /* register AMG type */
1506   if (!GAMGList) {
1507     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr);
1508     ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr);
1509   }
1510 
1511   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1512   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1513   pc->ops->setup          = PCSetUp_GAMG;
1514   pc->ops->reset          = PCReset_GAMG;
1515   pc->ops->destroy        = PCDestroy_GAMG;
1516 
1517   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1518                                             "PCGAMGSetProcEqLim_C",
1519                                             "PCGAMGSetProcEqLim_GAMG",
1520                                             PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1521   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1522                                             "PCGAMGSetCoarseEqLim_C",
1523                                             "PCGAMGSetCoarseEqLim_GAMG",
1524                                             PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1525   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1526                                             "PCGAMGSetRepartitioning_C",
1527                                             "PCGAMGSetRepartitioning_GAMG",
1528                                             PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr);
1529   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1530                                             "PCGAMGSetReuseProl_C",
1531                                             "PCGAMGSetReuseProl_GAMG",
1532                                             PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr);
1533   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1534                                             "PCGAMGSetUseASMAggs_C",
1535                                             "PCGAMGSetUseASMAggs_GAMG",
1536                                             PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr);
1537   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1538                                             "PCGAMGSetThreshold_C",
1539                                             "PCGAMGSetThreshold_GAMG",
1540                                             PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1541   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
1542                                             "PCGAMGSetType_C",
1543                                             "PCGAMGSetType_GAMG",
1544                                             PCGAMGSetType_GAMG);CHKERRQ(ierr);
1545   pc_gamg->repart = PETSC_FALSE;
1546   pc_gamg->reuse_prol = PETSC_TRUE;
1547   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1548   pc_gamg->min_eq_proc = 100;
1549   pc_gamg->coarse_eq_limit = 800;
1550   pc_gamg->threshold = 0.001;
1551   pc_gamg->Nlevels = GAMG_MAXLEVELS;
1552   pc_gamg->verbose = 0;
1553   pc_gamg->emax_id = -1;
1554   pc_gamg->eigtarget[0] = 0.05;
1555   pc_gamg->eigtarget[1] = 1.05;
1556 
1557   /* private events */
1558 #if defined PETSC_GAMG_USE_LOG
1559   if (count++ == 0) {
1560     ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
1561     ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
1562     /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1563     /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1564     /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1565     ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
1566     ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
1567     ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1568     ierr = PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
1569     ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
1570     ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
1571     ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
1572     ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
1573     ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
1574     ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
1575     ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
1576     ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
1577 
1578     /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1579     /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1580     /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1581     /* create timer stages */
1582 #if defined GAMG_STAGES
1583     {
1584       char     str[32];
1585       PetscInt lidx;
1586       sprintf(str,"MG Level %d (finest)",0);
1587       ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
1588       for (lidx=1;lidx<9;lidx++) {
1589         sprintf(str,"MG Level %d",lidx);
1590         ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
1591       }
1592     }
1593 #endif
1594   }
1595 #endif
1596   /* general events */
1597 #if defined PETSC_USE_LOG
1598   ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr);
1599   ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr);
1600   ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1601   ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1602   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1603   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1604   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr);
1605   ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr);
1606 #endif
1607 
1608   /* instantiate derived type */
1609   ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr);
1610   {
1611     char tname[256] = GAMGAGG;
1612     ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), PETSC_NULL);CHKERRQ(ierr);
1613     ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr);
1614   }
1615   ierr = PetscOptionsTail();CHKERRQ(ierr);
1616   PetscFunctionReturn(0);
1617 }
1618 EXTERN_C_END
1619