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