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