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