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