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