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