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