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