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