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