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