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