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