xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 4cdfd227e2bdf8feb38b6d90cfcac175440fcd00)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */
7f96513f1SMatthew G Knepley 
80cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
90cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
10b4fbaa2aSMark F. Adams #endif
110cbbd2e1SMark F. Adams 
120cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
13fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
150cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
19fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
200cbbd2e1SMark F. Adams #endif
210cbbd2e1SMark F. Adams 
22b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
230cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
24c1eae691SMark Adams static PetscLogStage gamg_stages[PETSC_GAMG_MAXLEVELS];
25b4fbaa2aSMark F. Adams #endif
26f96513f1SMatthew G Knepley 
27140e18c1SBarry Smith static PetscFunctionList GAMGList = 0;
283e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
299d5b6da9SMark F. Adams 
30d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
31d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
32d3d6bff4SMark F. Adams {
33d3d6bff4SMark F. Adams   PetscErrorCode ierr;
34d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
35d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
36d3d6bff4SMark F. Adams 
37d3d6bff4SMark F. Adams   PetscFunctionBegin;
3822a233eaSStefano Zampini   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
391c1aac46SBarry Smith   pc_gamg->data_sz = 0;
40878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
41a2f3521dSMark F. Adams   PetscFunctionReturn(0);
42a2f3521dSMark F. Adams }
43a2f3521dSMark F. Adams 
445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
455b89ad90SMark F. Adams /*
46c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
47a147abb0SMark F. Adams      of active processors.
485b89ad90SMark F. Adams 
495b89ad90SMark F. Adams    Input Parameter:
50a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
51a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
529d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
53c5bfad50SMark F. Adams    . cr_bs - coarse block size
543530afc2SMark F. Adams    In/Output Parameter:
55a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
56afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5711e60469SMark F. Adams    Output Parameter:
583530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
595b89ad90SMark F. Adams */
605cb416c2SMark F. Adams 
61171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
625b89ad90SMark F. Adams {
63a2f3521dSMark F. Adams   PetscErrorCode  ierr;
649d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
65486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
66a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
673b4367a7SBarry Smith   MPI_Comm        comm;
68c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
693ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
705b89ad90SMark F. Adams 
715b89ad90SMark F. Adams   PetscFunctionBegin;
723b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
733b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
743b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
75c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
769d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
77038e3b61SMark F. Adams 
78ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
79ce7c7f2fSMark Adams 
803ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
810298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
823ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
833ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
8473911c69SBarry Smith   } else {
853ae0bb68SMark Adams     PetscInt  bs;
863ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
873ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
883ae0bb68SMark Adams   }
89c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
90171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
91171cca9aSMark Adams   else {
92472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
930298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
94a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
957f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
96c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
97a2f3521dSMark F. Adams   }
98f852f58cSMark F. Adams 
992e3501ffSMark Adams   if (new_size==nactive) {
100ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
101ce7c7f2fSMark Adams     if (new_size < size) {
102ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
10331cb4603SMark Adams       ierr = PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr);
104ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
105ce7c7f2fSMark Adams         ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
106ce7c7f2fSMark Adams         ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
107ce7c7f2fSMark Adams       }
108ce7c7f2fSMark Adams     }
109ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1102e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
111192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
112885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
11371959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
11471959b99SBarry Smith     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);
115ce7c7f2fSMark Adams     /* get new_size and rfactor */
116ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
117ce7c7f2fSMark Adams       /* find factor */
118ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
119ce7c7f2fSMark Adams       else {
120ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
121ce7c7f2fSMark Adams         jj = -1;
122ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
123ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
124ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
125ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
126ce7c7f2fSMark Adams             if (fact > best_fact) {
127ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
128ce7c7f2fSMark Adams             }
129ce7c7f2fSMark Adams           }
130ce7c7f2fSMark Adams         }
131ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
132ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
133ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
134ce7c7f2fSMark Adams         else expand_factor = rfactor;
135ce7c7f2fSMark Adams       }
136ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
137*4cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
138*4cdfd227SMark         *a_Amat_crs = Cmat;
13931cb4603SMark Adams         ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
140ce7c7f2fSMark Adams         PetscFunctionReturn(0);
141ce7c7f2fSMark Adams       }
142ce7c7f2fSMark Adams     }
143*4cdfd227SMark #if defined PETSC_GAMG_USE_LOG
144*4cdfd227SMark     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
145*4cdfd227SMark #endif
146a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
147785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1482e3501ffSMark Adams     if (pc_gamg->repart) {
149a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1505a9b9e01SMark F. Adams       Mat      adj;
15131cb4603SMark Adams       ierr = PetscInfo4(pc,"Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");CHKERRQ(ierr);
152a2f3521dSMark F. Adams       /* get 'adj' */
153c5bfad50SMark F. Adams       if (cr_bs == 1) {
154038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
155806fa848SBarry Smith       } else {
156a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
157eb07cef2SMark F. Adams         Mat               tMat;
158a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
159b4fbaa2aSMark F. Adams         const PetscScalar *vals;
160b4fbaa2aSMark F. Adams         const PetscInt    *idx;
161a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
16239d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
163d9558ea9SBarry Smith         MatType           mtype;
164b4fbaa2aSMark F. Adams 
165e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
166a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
167a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
168c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
16958471d46SMark F. Adams           ierr      = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
170c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
171c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
17258471d46SMark F. Adams           ierr      = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr);
1733ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1743ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
17558471d46SMark F. Adams         }
1766876a03eSMark F. Adams 
177d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1783b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1793ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
180d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
181a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
182a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
183e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
184eb07cef2SMark F. Adams 
185a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
186c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
18722063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
188eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
189c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
190eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
191eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
192eb07cef2SMark F. Adams           }
19322063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
194eb07cef2SMark F. Adams         }
195eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
196eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197eb07cef2SMark F. Adams 
198b4fbaa2aSMark F. Adams         if (llev++ == -1) {
199b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2008caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2013b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
202b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2033bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
204b4fbaa2aSMark F. Adams         }
205eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
206eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
207a2f3521dSMark F. Adams       } /* create 'adj' */
208f150b916SMark F. Adams 
209a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2105a9b9e01SMark F. Adams         char            prefix[256];
2115a9b9e01SMark F. Adams         const char      *pcpre;
212b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
213b4fbaa2aSMark F. Adams         MatPartitioning mpart;
214a4b7d37bSMark F. Adams         IS              proc_is;
2152f03bc48SMark F. Adams 
2163b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2175ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2189d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2198caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
22059a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
22111e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
222c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
223a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
22411e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2255a9b9e01SMark F. Adams 
2265ef31b24SMark F. Adams         /* collect IS info */
227785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
228a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
229a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
230c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
231ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
232eb07cef2SMark F. Adams           }
2335ef31b24SMark F. Adams         }
234a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
235a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2365ef31b24SMark F. Adams       }
2375ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2385a9b9e01SMark F. Adams 
2393b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2408263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
24131cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
242ce7c7f2fSMark Adams       PetscInt targetPE;
243*4cdfd227SMark       if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
244302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
245ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2463b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
247a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
248e33ef3b1SMark F. Adams 
24911e60469SMark F. Adams     /*
250a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
25111e60469SMark F. Adams     */
252a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2537700e67bSMark Adams     is_eq_num_prim = is_eq_num;
25411e60469SMark F. Adams     /*
255a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
25611e60469SMark F. Adams     */
257c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
258c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
259a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
260ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
261a2f3521dSMark F. Adams 
262a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
263885364a3SMark Adams     /* 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 */
264885364a3SMark Adams     {
265885364a3SMark Adams       Vec            src_crd, dest_crd;
266885364a3SMark Adams       const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
267885364a3SMark Adams       VecScatter     vecscat;
268885364a3SMark Adams       PetscScalar    *array;
269885364a3SMark Adams       IS isscat;
270a2f3521dSMark F. Adams       /* move data (for primal equations only) */
27122063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
2723b4367a7SBarry Smith       ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2733ae0bb68SMark Adams       ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
274c0dedaeaSBarry Smith       ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
27511e60469SMark F. Adams       /*
2769d5b6da9SMark F. Adams 	There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
277c5bfad50SMark F. Adams 	a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
27811e60469SMark F. Adams       */
279854ce69bSBarry Smith       ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
280a2f3521dSMark F. Adams       ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2813ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
282c5bfad50SMark F. Adams 	PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
283a2f3521dSMark F. Adams 	for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
28411e60469SMark F. Adams       }
285a2f3521dSMark F. Adams       ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2863ae0bb68SMark Adams       ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
28792a756f0SMark F. Adams       ierr = PetscFree(tidx);CHKERRQ(ierr);
28811e60469SMark F. Adams       /*
28911e60469SMark F. Adams 	Create a vector to contain the original vertex information for each element
29011e60469SMark F. Adams       */
2913ae0bb68SMark Adams       ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2929d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
2933ae0bb68SMark Adams 	const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
2943ae0bb68SMark Adams 	for (ii=0; ii<ncrs; ii++) {
2959d5b6da9SMark F. Adams 	  for (kk=0; kk<ndata_rows; kk++) {
296a2f3521dSMark F. Adams 	    PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
297c8b0795cSMark F. Adams 	    PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
298676e1743SMark F. Adams 	    ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
299d3d6bff4SMark F. Adams 	  }
300038e3b61SMark F. Adams 	}
301eb07cef2SMark F. Adams       }
302eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
303eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
30411e60469SMark F. Adams       /*
30511e60469SMark F. Adams 	Scatter the element vertex information (still in the original vertex ordering)
30611e60469SMark F. Adams 	to the correct processor
30711e60469SMark F. Adams       */
3089448b7f1SJunchao Zhang       ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
30911e60469SMark F. Adams       ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31011e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31111e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31211e60469SMark F. Adams       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
31311e60469SMark F. Adams       ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
31411e60469SMark F. Adams       /*
31511e60469SMark F. Adams 	Put the element vertex data into a new allocation of the gdata->ele
31611e60469SMark F. Adams       */
317c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
318578f55a3SPeter Brune       ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3192fa5cd67SKarl Rupp 
3203ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3213ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3222fa5cd67SKarl Rupp 
32311e60469SMark F. Adams       ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3249d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3253ae0bb68SMark Adams 	for (ii=0; ii<ncrs_new; ii++) {
3269d5b6da9SMark F. Adams 	  for (kk=0; kk<ndata_rows; kk++) {
327a2f3521dSMark F. Adams 	    PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
328c8b0795cSMark F. Adams 	    pc_gamg->data[ix] = PetscRealPart(array[jx]);
329d3d6bff4SMark F. Adams 	  }
330038e3b61SMark F. Adams 	}
331038e3b61SMark F. Adams       }
33211e60469SMark F. Adams       ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
33311e60469SMark F. Adams       ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
334885364a3SMark Adams     }
335a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3370cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
338ed3f9983SMark F. Adams #endif
33911e60469SMark F. Adams     /*
3407dae84e0SHong Zhang       Invert for MatCreateSubMatrix
34111e60469SMark F. Adams     */
342a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
343a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
344c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
345a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
346a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
347a2f3521dSMark F. Adams     }
3483cb8563fSToby Isaac     if (Pcolumnperm) {
3493cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3503cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3513cb8563fSToby Isaac     }
352a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3530cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3540cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3550cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
356ed3f9983SMark F. Adams #endif
357a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
358a2f3521dSMark F. Adams     {
359a2f3521dSMark F. Adams       Mat mat;
3607dae84e0SHong Zhang       ierr        = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
361a2f3521dSMark F. Adams       *a_Amat_crs = mat;
362a2f3521dSMark F. Adams     }
363038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
364a2f3521dSMark F. Adams 
3650cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3660cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
367ed3f9983SMark F. Adams #endif
36811e60469SMark F. Adams     /* prolongator */
36911e60469SMark F. Adams     {
37011e60469SMark F. Adams       IS       findices;
371a2f3521dSMark F. Adams       PetscInt Istart,Iend;
372a2f3521dSMark F. Adams       Mat      Pnew;
37362294041SBarry Smith 
374a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3750cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3760cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
377ed3f9983SMark F. Adams #endif
3783b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
379c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
3807dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38111e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
382c5bfad50SMark F. Adams 
3830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3840cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
385ed3f9983SMark F. Adams #endif
3863530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
387a2f3521dSMark F. Adams 
388a2f3521dSMark F. Adams       /* output - repartitioned */
389a2f3521dSMark F. Adams       *a_P_inout = Pnew;
390e33ef3b1SMark F. Adams     }
391a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
3925b89ad90SMark F. Adams 
393c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
394ce7c7f2fSMark Adams 
395ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
396ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
397ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
3988bca76a6SMark Adams       static PetscInt llev = 2;
39939d09545SMark Adams       ierr = PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr);
400ce7c7f2fSMark Adams #endif
401ce7c7f2fSMark Adams       ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
402ce7c7f2fSMark Adams       ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
403ce7c7f2fSMark Adams       if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
404ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
405ce7c7f2fSMark Adams         PetscMPIInt size;
406ce7c7f2fSMark Adams         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
407ce7c7f2fSMark Adams         if (size > 1) {
408ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
409ce7c7f2fSMark Adams           ierr = VecPinToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr);
410ce7c7f2fSMark Adams           ierr = VecPinToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr);
411ce7c7f2fSMark Adams         }
412ce7c7f2fSMark Adams       }
413ce7c7f2fSMark Adams     }
414*4cdfd227SMark #if defined PETSC_GAMG_USE_LOG
415*4cdfd227SMark     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
416*4cdfd227SMark #endif
417a2f3521dSMark F. Adams   }
4185b89ad90SMark F. Adams   PetscFunctionReturn(0);
4195b89ad90SMark F. Adams }
4205b89ad90SMark F. Adams 
4215b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4225b89ad90SMark F. Adams /*
4235b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4245b89ad90SMark F. Adams                     by setting data structures and options.
4255b89ad90SMark F. Adams 
4265b89ad90SMark F. Adams    Input Parameter:
4275b89ad90SMark F. Adams .  pc - the preconditioner context
4285b89ad90SMark F. Adams 
4295b89ad90SMark F. Adams */
4309d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4315b89ad90SMark F. Adams {
4325b89ad90SMark F. Adams   PetscErrorCode ierr;
4339d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4345b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4352adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
436c1eae691SMark Adams   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_GAMG_MAXLEVELS];
4373b4367a7SBarry Smith   MPI_Comm       comm;
438c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
439c1eae691SMark Adams   Mat            Aarr[PETSC_GAMG_MAXLEVELS],Parr[PETSC_GAMG_MAXLEVELS];
440c1eae691SMark Adams   IS             *ASMLocalIDsArr[PETSC_GAMG_MAXLEVELS];
441a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
442569f4572SMark Adams   MatInfo        info;
443171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4445ef31b24SMark F. Adams 
4455b89ad90SMark F. Adams   PetscFunctionBegin;
4463b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4473b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4483b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
449dfd5c07aSMark F. Adams 
45084d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4511c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
452878e152fSMark F. Adams       /* reset everything */
453878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
454878e152fSMark F. Adams       pc->setupcalled = 0;
455806fa848SBarry Smith     } else {
45684d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
45703a628feSMark F. Adams       /* just do Galerkin grids */
45858471d46SMark F. Adams       Mat          B,dA,dB;
45958471d46SMark F. Adams 
46071959b99SBarry Smith       if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4619d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
46258471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
46323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
46458471d46SMark F. Adams         /* (re)set to get dirty flag */
46523ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
46658471d46SMark F. Adams 
4672fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
468ef3f0257SMark Adams           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
469ef3f0257SMark Adams           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
470ef3f0257SMark Adams             ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
471ef3f0257SMark Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr);
472084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
47303a628feSMark F. Adams             mglevels[level]->A = B;
474806fa848SBarry Smith           } else {
475ef3f0257SMark Adams             ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
47623ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
47758471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
47803a628feSMark F. Adams           }
47923ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
48058471d46SMark F. Adams           dB   = B;
48158471d46SMark F. Adams         }
4825f8cf99dSMark F. Adams       }
483d5280255SMark F. Adams 
484d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
48558471d46SMark F. Adams       PetscFunctionReturn(0);
486eb07cef2SMark F. Adams     }
487878e152fSMark F. Adams   }
488f6536408SMark F. Adams 
489878e152fSMark F. Adams   if (!pc_gamg->data) {
490878e152fSMark F. Adams     if (pc_gamg->orig_data) {
491878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4920298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
4932fa5cd67SKarl Rupp 
494878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
495878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
496878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
4972fa5cd67SKarl Rupp 
498785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
499878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
500806fa848SBarry Smith     } else {
5011ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5027700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5039d5b6da9SMark F. Adams     }
504878e152fSMark F. Adams   }
505878e152fSMark F. Adams 
506878e152fSMark F. Adams   /* cache original data for reuse */
5071c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
508785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
509878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
510878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
511878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
512878e152fSMark F. Adams   }
513038e3b61SMark F. Adams 
514302f38e8SMark F. Adams   /* get basic dims */
515302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
516171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
51784d3f75bSMark F. Adams 
518569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
519569f4572SMark Adams   nnz0   = info.nz_used;
520569f4572SMark Adams   nnztot = info.nz_used;
52162294041SBarry Smith   ierr = PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr);
522569f4572SMark Adams 
523a2f3521dSMark F. Adams   /* Get A_i and R_i */
52462294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5259ab59c8bSMark Adams     pc_gamg->current_level = level;
5265b89ad90SMark F. Adams     level1 = level + 1;
5270cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5280cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
529a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
530a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
531b4fbaa2aSMark F. Adams #endif
532a2f3521dSMark F. Adams #endif
533c8b0795cSMark F. Adams     { /* construct prolongator */
534725b86d8SJed Brown       Mat              Gmat;
5350cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5367700e67bSMark Adams       Mat              Prol11;
537c8b0795cSMark F. Adams 
5387700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5391ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5407700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
541c8b0795cSMark F. Adams 
542a2f3521dSMark F. Adams       /* could have failed to create new level */
543a2f3521dSMark F. Adams       if (Prol11) {
5449d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5450298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
546a2f3521dSMark F. Adams 
547fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
548c8b0795cSMark F. Adams           /* smooth */
549fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
550c8b0795cSMark F. Adams         }
551c8b0795cSMark F. Adams 
5520c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
5531b18a24aSMark Adams           PetscInt bs;
5541b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5550a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
556ffc955d6SMark F. Adams         }
557ffc955d6SMark F. Adams 
5584bde40a0SMark Adams         Parr[level1] = Prol11;
5594bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
5604bde40a0SMark Adams 
561a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
56241b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
563a2f3521dSMark F. Adams     } /* construct prolongator scope */
5640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5650cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
566c8b0795cSMark F. Adams #endif
5677f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
568171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
569569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
57062294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
571a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
572a90e85d9SMark Adams #endif
573c8b0795cSMark F. Adams       break;
574c8b0795cSMark F. Adams     }
5750cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5760cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
577b4fbaa2aSMark F. Adams #endif
578171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
579171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
580171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
5810e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
582171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
583a2f3521dSMark F. Adams 
5840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5850cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
586b4fbaa2aSMark F. Adams #endif
587171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
588569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
589569f4572SMark Adams     nnztot += info.nz_used;
5901d5b2942SMark Adams     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);
591569f4572SMark Adams 
5920cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
593b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
594b4fbaa2aSMark F. Adams #endif
595a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
5969ab59c8bSMark Adams     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
5979ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
598a90e85d9SMark Adams       level++;
599a90e85d9SMark Adams       break;
600a90e85d9SMark Adams     }
601c8b0795cSMark F. Adams   } /* levels */
602c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
603c8b0795cSMark F. Adams 
604569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
6059d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6065b89ad90SMark F. Adams   fine_level       = level;
6070298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6085b89ad90SMark F. Adams 
60962294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
610d5280255SMark F. Adams     /* set default smoothers & set operators */
61162294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
612ffc955d6SMark F. Adams       KSP smoother;
613ffc955d6SMark F. Adams       PC  subpc;
614a2f3521dSMark F. Adams 
6159d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
616f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
617ffc955d6SMark F. Adams 
618a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
619a2f3521dSMark F. Adams       /* set ops */
62023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
621a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
622a2f3521dSMark F. Adams 
623a2f3521dSMark F. Adams       /* set defaults */
6246c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
625a2f3521dSMark F. Adams 
6260c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6270c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6282d3561bbSSatish Balay         PetscInt sz;
6297a28f3e5SMark Adams         IS       *iss;
630a2f3521dSMark F. Adams 
6312d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6327a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
6330c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6340a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6350c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6367f66b68fSMark Adams         if (!sz) {
637ffc955d6SMark F. Adams           IS       is;
6380a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6397a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
640a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
641806fa848SBarry Smith         } else {
642a94c3b12SMark F. Adams           PetscInt kk;
6437a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
644a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6457a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
646a94c3b12SMark F. Adams           }
6477a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
648ffc955d6SMark F. Adams         }
6490298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
650ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
651806fa848SBarry Smith       } else {
652890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
653ffc955d6SMark F. Adams       }
654d5280255SMark F. Adams     }
655d5280255SMark F. Adams     {
656d5280255SMark F. Adams       /* coarse grid */
657d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
658d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
659d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
66023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
661cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
662d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
663d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
664d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
665d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
66671959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
66771959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
668d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
669d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
6709dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
6712fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
67208e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
6735b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
6745b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
6755b42dca8SJed Brown          * view every subdomain as though they were different. */
6765b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
677d5280255SMark F. Adams       }
678cf8ae1d3SMark Adams     }
679d5280255SMark F. Adams 
680d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
681d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
682e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
683d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
68469aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
685d5280255SMark F. Adams 
686d5280255SMark F. Adams     /* clean up */
687d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
688587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
689587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
6905b89ad90SMark F. Adams     }
6910cbbd2e1SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
692806fa848SBarry Smith   } else {
6935f8cf99dSMark F. Adams     KSP smoother;
694302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
6959d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
69623ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
6975f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
6989d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6995f8cf99dSMark F. Adams   }
7005b89ad90SMark F. Adams   PetscFunctionReturn(0);
7015b89ad90SMark F. Adams }
7025b89ad90SMark F. Adams 
703eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7045b89ad90SMark F. Adams /*
7055b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7065b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7075b89ad90SMark F. Adams 
7085b89ad90SMark F. Adams    Input Parameter:
7095b89ad90SMark F. Adams .  pc - the preconditioner context
7105b89ad90SMark F. Adams 
7115b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7125b89ad90SMark F. Adams */
7135b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7145b89ad90SMark F. Adams {
7155b89ad90SMark F. Adams   PetscErrorCode ierr;
7165b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7175b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7185b89ad90SMark F. Adams 
7195b89ad90SMark F. Adams   PetscFunctionBegin;
7205b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7219b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7229b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7239b8ffb57SJed Brown   }
7241ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7251ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7265b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7275b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7285b89ad90SMark F. Adams   PetscFunctionReturn(0);
7295b89ad90SMark F. Adams }
7305b89ad90SMark F. Adams 
731676e1743SMark F. Adams /*@
732cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
733676e1743SMark F. Adams 
7341cc46a46SBarry Smith    Logically Collective on PC
735676e1743SMark F. Adams 
736676e1743SMark F. Adams    Input Parameters:
7371cc46a46SBarry Smith +  pc - the preconditioner context
7381cc46a46SBarry Smith -  n - the number of equations
739676e1743SMark F. Adams 
740676e1743SMark F. Adams 
741676e1743SMark F. Adams    Options Database Key:
7421cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
743676e1743SMark F. Adams 
74495452b02SPatrick Sanan    Notes:
74595452b02SPatrick Sanan     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
746cab9ed1eSBarry Smith           that has degrees of freedom
747cab9ed1eSBarry Smith 
748676e1743SMark F. Adams    Level: intermediate
749676e1743SMark F. Adams 
7501c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
751676e1743SMark F. Adams @*/
752676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
753676e1743SMark F. Adams {
754676e1743SMark F. Adams   PetscErrorCode ierr;
755676e1743SMark F. Adams 
756676e1743SMark F. Adams   PetscFunctionBegin;
757676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
758676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
759676e1743SMark F. Adams   PetscFunctionReturn(0);
760676e1743SMark F. Adams }
761676e1743SMark F. Adams 
7621e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
763676e1743SMark F. Adams {
764c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
765c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
766676e1743SMark F. Adams 
767676e1743SMark F. Adams   PetscFunctionBegin;
7689d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
769676e1743SMark F. Adams   PetscFunctionReturn(0);
770676e1743SMark F. Adams }
771676e1743SMark F. Adams 
772389730f3SMark F. Adams /*@
773cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
774389730f3SMark F. Adams 
775389730f3SMark F. Adams  Collective on PC
776389730f3SMark F. Adams 
777389730f3SMark F. Adams    Input Parameters:
7781cc46a46SBarry Smith +  pc - the preconditioner context
7791cc46a46SBarry Smith -  n - maximum number of equations to aim for
780389730f3SMark F. Adams 
781389730f3SMark F. Adams    Options Database Key:
7821cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
783389730f3SMark F. Adams 
78474329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
78574329af1SBarry Smith      has less than 1000 unknowns.
78674329af1SBarry Smith 
787389730f3SMark F. Adams    Level: intermediate
788389730f3SMark F. Adams 
7891c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
790389730f3SMark F. Adams @*/
791389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
792389730f3SMark F. Adams {
793389730f3SMark F. Adams   PetscErrorCode ierr;
794389730f3SMark F. Adams 
795389730f3SMark F. Adams   PetscFunctionBegin;
796389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
797389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
798389730f3SMark F. Adams   PetscFunctionReturn(0);
799389730f3SMark F. Adams }
800389730f3SMark F. Adams 
8011e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
802389730f3SMark F. Adams {
803389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
804389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
805389730f3SMark F. Adams 
806389730f3SMark F. Adams   PetscFunctionBegin;
8079d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
808389730f3SMark F. Adams   PetscFunctionReturn(0);
809389730f3SMark F. Adams }
810389730f3SMark F. Adams 
811676e1743SMark F. Adams /*@
812cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
813676e1743SMark F. Adams 
814676e1743SMark F. Adams    Collective on PC
815676e1743SMark F. Adams 
816676e1743SMark F. Adams    Input Parameters:
8171cc46a46SBarry Smith +  pc - the preconditioner context
8181cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
819676e1743SMark F. Adams 
820676e1743SMark F. Adams    Options Database Key:
8211cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
822676e1743SMark F. Adams 
82395452b02SPatrick Sanan    Notes:
82495452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
825cab9ed1eSBarry Smith 
826676e1743SMark F. Adams    Level: intermediate
827676e1743SMark F. Adams 
828676e1743SMark F. Adams .seealso: ()
829676e1743SMark F. Adams @*/
830cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
831676e1743SMark F. Adams {
832676e1743SMark F. Adams   PetscErrorCode ierr;
833676e1743SMark F. Adams 
834676e1743SMark F. Adams   PetscFunctionBegin;
835676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
836cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
837676e1743SMark F. Adams   PetscFunctionReturn(0);
838676e1743SMark F. Adams }
839676e1743SMark F. Adams 
840cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
841676e1743SMark F. Adams {
842c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
843c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
844676e1743SMark F. Adams 
845676e1743SMark F. Adams   PetscFunctionBegin;
8469d5b6da9SMark F. Adams   pc_gamg->repart = n;
847676e1743SMark F. Adams   PetscFunctionReturn(0);
848676e1743SMark F. Adams }
849676e1743SMark F. Adams 
850dfd5c07aSMark F. Adams /*@
851cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
852dfd5c07aSMark F. Adams 
853dfd5c07aSMark F. Adams    Collective on PC
854dfd5c07aSMark F. Adams 
855dfd5c07aSMark F. Adams    Input Parameters:
8561cc46a46SBarry Smith +  pc - the preconditioner context
8571cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
858dfd5c07aSMark F. Adams 
859dfd5c07aSMark F. Adams    Options Database Key:
8601cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
861dfd5c07aSMark F. Adams 
862dfd5c07aSMark F. Adams    Level: intermediate
863dfd5c07aSMark F. Adams 
86495452b02SPatrick Sanan    Notes:
86595452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
866cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
867cab9ed1eSBarry Smith 
868dfd5c07aSMark F. Adams .seealso: ()
869dfd5c07aSMark F. Adams @*/
8701cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
871dfd5c07aSMark F. Adams {
872dfd5c07aSMark F. Adams   PetscErrorCode ierr;
873dfd5c07aSMark F. Adams 
874dfd5c07aSMark F. Adams   PetscFunctionBegin;
875dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8761cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
877dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
878dfd5c07aSMark F. Adams }
879dfd5c07aSMark F. Adams 
8801cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
881dfd5c07aSMark F. Adams {
882dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
883dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
884dfd5c07aSMark F. Adams 
885dfd5c07aSMark F. Adams   PetscFunctionBegin;
886dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
887dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
888dfd5c07aSMark F. Adams }
889dfd5c07aSMark F. Adams 
890ffc955d6SMark F. Adams /*@
891cab9ed1eSBarry Smith    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
892ffc955d6SMark F. Adams 
893ffc955d6SMark F. Adams    Collective on PC
894ffc955d6SMark F. Adams 
895ffc955d6SMark F. Adams    Input Parameters:
896cab9ed1eSBarry Smith +  pc - the preconditioner context
897cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
898ffc955d6SMark F. Adams 
899ffc955d6SMark F. Adams    Options Database Key:
900cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
901ffc955d6SMark F. Adams 
902ffc955d6SMark F. Adams    Level: intermediate
903ffc955d6SMark F. Adams 
904ffc955d6SMark F. Adams .seealso: ()
905ffc955d6SMark F. Adams @*/
906cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
907ffc955d6SMark F. Adams {
908ffc955d6SMark F. Adams   PetscErrorCode ierr;
909ffc955d6SMark F. Adams 
910ffc955d6SMark F. Adams   PetscFunctionBegin;
911ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
912cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
913ffc955d6SMark F. Adams   PetscFunctionReturn(0);
914ffc955d6SMark F. Adams }
915ffc955d6SMark F. Adams 
916cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
917ffc955d6SMark F. Adams {
918ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
919ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
920ffc955d6SMark F. Adams 
921ffc955d6SMark F. Adams   PetscFunctionBegin;
922cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
923ffc955d6SMark F. Adams   PetscFunctionReturn(0);
924ffc955d6SMark F. Adams }
925ffc955d6SMark F. Adams 
926171cca9aSMark Adams /*@
927cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
928171cca9aSMark Adams 
929171cca9aSMark Adams    Collective on PC
930171cca9aSMark Adams 
931171cca9aSMark Adams    Input Parameters:
932171cca9aSMark Adams +  pc - the preconditioner context
933cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
934171cca9aSMark Adams 
935171cca9aSMark Adams    Options Database Key:
936cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
937171cca9aSMark Adams 
938171cca9aSMark Adams    Level: intermediate
939171cca9aSMark Adams 
94039d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
941171cca9aSMark Adams @*/
942171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
943171cca9aSMark Adams {
944171cca9aSMark Adams   PetscErrorCode ierr;
945171cca9aSMark Adams 
946171cca9aSMark Adams   PetscFunctionBegin;
947171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
948171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
949171cca9aSMark Adams   PetscFunctionReturn(0);
950171cca9aSMark Adams }
951171cca9aSMark Adams 
952171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
953171cca9aSMark Adams {
954171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
955171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
956171cca9aSMark Adams 
957171cca9aSMark Adams   PetscFunctionBegin;
958171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
959ffc955d6SMark F. Adams   PetscFunctionReturn(0);
960ffc955d6SMark F. Adams }
961ffc955d6SMark F. Adams 
9624ef23d27SMark F. Adams /*@
963ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
964ce7c7f2fSMark Adams 
965ce7c7f2fSMark Adams    Collective on PC
966ce7c7f2fSMark Adams 
967ce7c7f2fSMark Adams    Input Parameters:
968ce7c7f2fSMark Adams +  pc - the preconditioner context
969ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
970ce7c7f2fSMark Adams 
971ce7c7f2fSMark Adams    Options Database Key:
972ce7c7f2fSMark Adams .  -pc_gamg_cpu_pin_coarse_grids
973ce7c7f2fSMark Adams 
974ce7c7f2fSMark Adams    Level: intermediate
975ce7c7f2fSMark Adams 
97639d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
977ce7c7f2fSMark Adams @*/
978ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
979ce7c7f2fSMark Adams {
980ce7c7f2fSMark Adams   PetscErrorCode ierr;
981ce7c7f2fSMark Adams 
982ce7c7f2fSMark Adams   PetscFunctionBegin;
983ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
984ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
985ce7c7f2fSMark Adams   PetscFunctionReturn(0);
986ce7c7f2fSMark Adams }
987ce7c7f2fSMark Adams 
988ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
989ce7c7f2fSMark Adams {
990ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
991ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
992ce7c7f2fSMark Adams 
993ce7c7f2fSMark Adams   PetscFunctionBegin;
994ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
995ce7c7f2fSMark Adams   PetscFunctionReturn(0);
996ce7c7f2fSMark Adams }
997ce7c7f2fSMark Adams 
998ce7c7f2fSMark Adams /*@
999ce7c7f2fSMark Adams    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1000ce7c7f2fSMark Adams 
1001ce7c7f2fSMark Adams    Collective on PC
1002ce7c7f2fSMark Adams 
1003ce7c7f2fSMark Adams    Input Parameters:
1004ce7c7f2fSMark Adams +  pc - the preconditioner context
1005ce7c7f2fSMark Adams -  flg - Layout type
1006ce7c7f2fSMark Adams 
1007ce7c7f2fSMark Adams    Options Database Key:
1008ce7c7f2fSMark Adams .  -pc_gamg_coarse_grid_layout_type
1009ce7c7f2fSMark Adams 
1010ce7c7f2fSMark Adams    Level: intermediate
1011ce7c7f2fSMark Adams 
101239d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1013ce7c7f2fSMark Adams @*/
1014ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1015ce7c7f2fSMark Adams {
1016ce7c7f2fSMark Adams   PetscErrorCode ierr;
1017ce7c7f2fSMark Adams 
1018ce7c7f2fSMark Adams   PetscFunctionBegin;
1019ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1020ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1021ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1022ce7c7f2fSMark Adams }
1023ce7c7f2fSMark Adams 
1024ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1025ce7c7f2fSMark Adams {
1026ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1027ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1028ce7c7f2fSMark Adams 
1029ce7c7f2fSMark Adams   PetscFunctionBegin;
1030ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1031ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1032ce7c7f2fSMark Adams }
1033ce7c7f2fSMark Adams 
1034ce7c7f2fSMark Adams /*@
10351cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
10364ef23d27SMark F. Adams 
10374ef23d27SMark F. Adams    Not collective on PC
10384ef23d27SMark F. Adams 
10394ef23d27SMark F. Adams    Input Parameters:
10401cc46a46SBarry Smith +  pc - the preconditioner
10411cc46a46SBarry Smith -  n - the maximum number of levels to use
10424ef23d27SMark F. Adams 
10434ef23d27SMark F. Adams    Options Database Key:
10444ef23d27SMark F. Adams .  -pc_mg_levels
10454ef23d27SMark F. Adams 
10464ef23d27SMark F. Adams    Level: intermediate
10474ef23d27SMark F. Adams 
10484ef23d27SMark F. Adams .seealso: ()
10494ef23d27SMark F. Adams @*/
10504ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
10514ef23d27SMark F. Adams {
10524ef23d27SMark F. Adams   PetscErrorCode ierr;
10534ef23d27SMark F. Adams 
10544ef23d27SMark F. Adams   PetscFunctionBegin;
10554ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10564ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
10574ef23d27SMark F. Adams   PetscFunctionReturn(0);
10584ef23d27SMark F. Adams }
10594ef23d27SMark F. Adams 
10601e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
10614ef23d27SMark F. Adams {
10624ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
10634ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
10644ef23d27SMark F. Adams 
10654ef23d27SMark F. Adams   PetscFunctionBegin;
10669d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
10674ef23d27SMark F. Adams   PetscFunctionReturn(0);
10684ef23d27SMark F. Adams }
10694ef23d27SMark F. Adams 
10703542efc5SMark F. Adams /*@
10713542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
10723542efc5SMark F. Adams 
10733542efc5SMark F. Adams    Not collective on PC
10743542efc5SMark F. Adams 
10753542efc5SMark F. Adams    Input Parameters:
10761cc46a46SBarry Smith +  pc - the preconditioner context
1077b001cb0fSBarry Smith -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
10783542efc5SMark F. Adams 
10793542efc5SMark F. Adams    Options Database Key:
10801cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
10813542efc5SMark F. Adams 
108295452b02SPatrick Sanan    Notes:
1083af3c827dSMark Adams     Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
1084af3c827dSMark Adams     Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.
1085cab9ed1eSBarry Smith 
10863542efc5SMark F. Adams    Level: intermediate
10873542efc5SMark F. Adams 
1088af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
10893542efc5SMark F. Adams @*/
1090c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
10913542efc5SMark F. Adams {
10923542efc5SMark F. Adams   PetscErrorCode ierr;
10933542efc5SMark F. Adams 
10943542efc5SMark F. Adams   PetscFunctionBegin;
10953542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1096c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
10973542efc5SMark F. Adams   PetscFunctionReturn(0);
10983542efc5SMark F. Adams }
10993542efc5SMark F. Adams 
1100c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
11013542efc5SMark F. Adams {
1102c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1103c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1104c1eae691SMark Adams   PetscInt i;
1105c1eae691SMark Adams   PetscFunctionBegin;
1106c1eae691SMark Adams   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1107c1eae691SMark Adams   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1108c1eae691SMark Adams   PetscFunctionReturn(0);
1109c1eae691SMark Adams }
1110c1eae691SMark Adams 
1111c1eae691SMark Adams /*@
1112c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1113c1eae691SMark Adams 
1114c1eae691SMark Adams    Not collective on PC
1115c1eae691SMark Adams 
1116c1eae691SMark Adams    Input Parameters:
1117c1eae691SMark Adams +  pc - the preconditioner context
1118c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1119c1eae691SMark Adams 
1120c1eae691SMark Adams    Options Database Key:
1121c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1122c1eae691SMark Adams 
1123c1eae691SMark Adams    Level: advanced
1124c1eae691SMark Adams 
1125c1eae691SMark Adams .seealso: ()
1126c1eae691SMark Adams @*/
1127c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1128c1eae691SMark Adams {
1129c1eae691SMark Adams   PetscErrorCode ierr;
11303542efc5SMark F. Adams 
11313542efc5SMark F. Adams   PetscFunctionBegin;
1132c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1133c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1134c1eae691SMark Adams   PetscFunctionReturn(0);
1135c1eae691SMark Adams }
1136c1eae691SMark Adams 
1137c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1138c1eae691SMark Adams {
1139c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1140c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1141c1eae691SMark Adams   PetscFunctionBegin;
1142c1eae691SMark Adams   pc_gamg->threshold_scale = v;
11433542efc5SMark F. Adams   PetscFunctionReturn(0);
11443542efc5SMark F. Adams }
11453542efc5SMark F. Adams 
1146e20c40e8SBarry Smith /*@C
1147c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1148676e1743SMark F. Adams 
1149676e1743SMark F. Adams    Collective on PC
1150676e1743SMark F. Adams 
1151676e1743SMark F. Adams    Input Parameters:
1152c60c7ad4SBarry Smith +  pc - the preconditioner context
1153c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1154676e1743SMark F. Adams 
1155676e1743SMark F. Adams    Options Database Key:
1156cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1157676e1743SMark F. Adams 
1158676e1743SMark F. Adams    Level: intermediate
1159676e1743SMark F. Adams 
1160cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1161676e1743SMark F. Adams @*/
116219fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1163676e1743SMark F. Adams {
1164676e1743SMark F. Adams   PetscErrorCode ierr;
1165676e1743SMark F. Adams 
1166676e1743SMark F. Adams   PetscFunctionBegin;
1167676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1168806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1169676e1743SMark F. Adams   PetscFunctionReturn(0);
1170676e1743SMark F. Adams }
1171676e1743SMark F. Adams 
1172e20c40e8SBarry Smith /*@C
1173c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1174c60c7ad4SBarry Smith 
1175c60c7ad4SBarry Smith    Collective on PC
1176c60c7ad4SBarry Smith 
1177c60c7ad4SBarry Smith    Input Parameter:
1178c60c7ad4SBarry Smith .  pc - the preconditioner context
1179c60c7ad4SBarry Smith 
1180c60c7ad4SBarry Smith    Output Parameter:
1181c60c7ad4SBarry Smith .  type - the type of algorithm used
1182c60c7ad4SBarry Smith 
1183c60c7ad4SBarry Smith    Level: intermediate
1184c60c7ad4SBarry Smith 
11851c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1186c60c7ad4SBarry Smith @*/
1187c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1188c60c7ad4SBarry Smith {
1189c60c7ad4SBarry Smith   PetscErrorCode ierr;
1190c60c7ad4SBarry Smith 
1191c60c7ad4SBarry Smith   PetscFunctionBegin;
1192c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1193c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1194c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1195c60c7ad4SBarry Smith }
1196c60c7ad4SBarry Smith 
1197c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1198c60c7ad4SBarry Smith {
1199c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1200c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1201c60c7ad4SBarry Smith 
1202c60c7ad4SBarry Smith   PetscFunctionBegin;
1203c60c7ad4SBarry Smith   *type = pc_gamg->type;
1204c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1205c60c7ad4SBarry Smith }
1206c60c7ad4SBarry Smith 
12071e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1208676e1743SMark F. Adams {
12099d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
12101ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
12111ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1212676e1743SMark F. Adams 
1213676e1743SMark F. Adams   PetscFunctionBegin;
1214c60c7ad4SBarry Smith   pc_gamg->type = type;
12151c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
12169d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
12171ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
12181ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
12191ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1220e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
12213ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
12223ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
12233ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
12243ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
12253ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
12263ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
12273ae0bb68SMark Adams     pc_gamg->data_sz = 0;
12281ab5ffc9SJed Brown   }
12291ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
12301ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
12319d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1232676e1743SMark F. Adams   PetscFunctionReturn(0);
1233676e1743SMark F. Adams }
1234676e1743SMark F. Adams 
12359d766c59SMark Adams /* -------------------------------------------------------------------------- */
12369d766c59SMark Adams /*
12379d766c59SMark Adams    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
12389d766c59SMark Adams 
12399d766c59SMark Adams    Input Parameter:
12409d766c59SMark Adams .  pc - the preconditioner context
12419d766c59SMark Adams 
12429d766c59SMark Adams    Output Parameter:
12439d766c59SMark Adams .  gc - grid complexity = sum_i(nnz_i) / nnz_0
12449d766c59SMark Adams 
12459d766c59SMark Adams    Level: advanced
12469d766c59SMark Adams */
12479d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
12489d766c59SMark Adams {
12499d766c59SMark Adams   PetscErrorCode ierr;
12509d766c59SMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
12519d766c59SMark Adams   PC_MG_Levels   **mglevels = mg->levels;
12529d766c59SMark Adams   PetscInt       lev, nnz0 = -1;
12539d766c59SMark Adams   MatInfo        info;
12549d766c59SMark Adams   PetscFunctionBegin;
12559d766c59SMark Adams   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
12569d766c59SMark Adams   for (lev=0, *gc=0; lev<mg->nlevels; lev++) {
125762a6e064SMark Adams     Mat dB;
125862a6e064SMark Adams     ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr);
125962a6e064SMark Adams     ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
12609d766c59SMark Adams     *gc += (PetscReal)info.nz_used;
12619d766c59SMark Adams     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
12629d766c59SMark Adams   }
12639d766c59SMark Adams   if (nnz0) *gc /= (PetscReal)nnz0;
12649d766c59SMark Adams   else *gc = 0;
12659d766c59SMark Adams   PetscFunctionReturn(0);
12669d766c59SMark Adams }
12679d766c59SMark Adams 
12685adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
12695adeb434SBarry Smith {
1270c1eae691SMark Adams   PetscErrorCode ierr,i;
12715adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
12725adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
127323b2d91dSMark Adams   PetscReal       gc=0;
12745adeb434SBarry Smith   PetscFunctionBegin;
12755adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1276459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1277c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1278c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1279c1eae691SMark Adams   }
1280459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1281459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1282cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1283cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1284cab9ed1eSBarry Smith   }
1285171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1286171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1287171cca9aSMark Adams   }
1288ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1289ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1290ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1291ce7c7f2fSMark Adams   }
1292ce7c7f2fSMark Adams #endif
1293ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1294ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1295ce7c7f2fSMark Adams   /* } else { */
1296ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1297ce7c7f2fSMark Adams   /* } */
12985adeb434SBarry Smith   if (pc_gamg->ops->view) {
12995adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
13005adeb434SBarry Smith   }
13019d766c59SMark Adams   ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr);
13029d766c59SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);CHKERRQ(ierr);
13035adeb434SBarry Smith   PetscFunctionReturn(0);
13045adeb434SBarry Smith }
13055adeb434SBarry Smith 
13064416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
13075b89ad90SMark F. Adams {
1308676e1743SMark F. Adams   PetscErrorCode ierr;
1309676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1310676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1311676e1743SMark F. Adams   PetscBool      flag;
13123b4367a7SBarry Smith   MPI_Comm       comm;
131314a9496bSBarry Smith   char           prefix[256];
1314c1eae691SMark Adams   PetscInt       i,n;
131514a9496bSBarry Smith   const char     *pcpre;
1316d605200eSMark Adams   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",0};
13175b89ad90SMark F. Adams   PetscFunctionBegin;
13183b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1319e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
1320676e1743SMark F. Adams   {
1321bd94a7aaSJed Brown     char tname[256];
13221a1c1e04SBarry Smith     ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1323bd94a7aaSJed Brown     if (flag) {
1324bd94a7aaSJed Brown       ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
13251ab5ffc9SJed Brown     }
1326cab9ed1eSBarry Smith     ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
13271cc46a46SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1328a303c832SJed Brown     ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr);
1329cf8ae1d3SMark Adams     ierr = PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);CHKERRQ(ierr);
1330ce7c7f2fSMark Adams     ierr = PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL);CHKERRQ(ierr);
1331a0095786SMark     ierr = PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL);CHKERRQ(ierr);
133294ae4db5SBarry Smith     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);
133394ae4db5SBarry Smith     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);
1334a303c832SJed Brown     ierr = PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);CHKERRQ(ierr);
1335c1eae691SMark Adams     n = PETSC_GAMG_MAXLEVELS;
1336c1eae691SMark Adams     ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
1337efd3c5ceSMark Adams     if (!flag || n < PETSC_GAMG_MAXLEVELS) {
1338efd3c5ceSMark Adams       if (!flag) n = 1;
1339c1eae691SMark Adams       i = n;
1340c1eae691SMark Adams       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1341c1eae691SMark Adams     }
134294ae4db5SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
1343b7cbab4eSMark Adams 
1344b7cbab4eSMark Adams     /* set options for subtype */
1345e55864a3SBarry Smith     if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
1346676e1743SMark F. Adams   }
134714a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
134814a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1349676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
13505b89ad90SMark F. Adams   PetscFunctionReturn(0);
13515b89ad90SMark F. Adams }
13525b89ad90SMark F. Adams 
13535b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
13545b89ad90SMark F. Adams /*MC
13551cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
13565b89ad90SMark F. Adams 
1357280d9858SJed Brown    Options Database Keys:
1358cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1359cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1360cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1361cab9ed1eSBarry Smith .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1362cab9ed1eSBarry Smith .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1363cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1364cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
13656008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1366c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1367cab9ed1eSBarry Smith 
1368cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1369cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1370cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1371cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1372cab9ed1eSBarry Smith 
1373db9745e2SBarry Smith    Multigrid options:
1374db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1375db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1376db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1377cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
13785b89ad90SMark F. Adams 
13791cc46a46SBarry Smith 
138095452b02SPatrick Sanan   Notes:
138195452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1382db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1383db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1384db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
13851cc46a46SBarry Smith 
13865b89ad90SMark F. Adams   Level: intermediate
1387280d9858SJed Brown 
13881cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1389171cca9aSMark Adams            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
13905b89ad90SMark F. Adams M*/
1391b2573a8aSBarry Smith 
13928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
13935b89ad90SMark F. Adams {
1394c1eae691SMark Adams   PetscErrorCode ierr,i;
13955b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
13965b89ad90SMark F. Adams   PC_MG          *mg;
13975b89ad90SMark F. Adams 
13985b89ad90SMark F. Adams   PetscFunctionBegin;
13991c1aac46SBarry Smith    /* register AMG type */
14001c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
14011c1aac46SBarry Smith 
14025b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
14031c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
14045b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
14055b89ad90SMark F. Adams 
14065b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1407b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
140869aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
14095b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
14105b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
14115b89ad90SMark F. Adams 
1412b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
14131ab5ffc9SJed Brown 
14149d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
14159d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
14169d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
14179d5b6da9SMark F. Adams   pc_gamg->data    = 0;
14185b89ad90SMark F. Adams 
14199d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
14205b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
14215b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
14225b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
14235b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
14245adeb434SBarry Smith   mg->view                = PCView_GAMG;
14255b89ad90SMark F. Adams 
142697d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
142797d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1428bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1429bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1430cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
14311cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1432cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1433171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1434ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1435ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1436bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1437c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1438bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1439c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1440bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
14419d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1442d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
14430c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1444171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1445a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1446a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1447038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
144825a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
1449c1eae691SMark Adams   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1450c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
1451c1eae691SMark Adams   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
14529ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1453c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14549d5b6da9SMark F. Adams 
1455bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1456bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
14575b89ad90SMark F. Adams   PetscFunctionReturn(0);
14585b89ad90SMark F. Adams }
14593e3471ccSMark Adams 
14603e3471ccSMark Adams /*@C
14613e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
14628a690491SBarry Smith     from PCInitializePackage().
14633e3471ccSMark Adams 
14643e3471ccSMark Adams  Level: developer
14653e3471ccSMark Adams 
14663e3471ccSMark Adams  .seealso: PetscInitialize()
14673e3471ccSMark Adams @*/
14683e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
14693e3471ccSMark Adams {
14703e3471ccSMark Adams   PetscErrorCode ierr;
14713e3471ccSMark Adams 
14723e3471ccSMark Adams   PetscFunctionBegin;
14733e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
14743e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
14753e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
14763e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
14778e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
14783e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1479c1c463dbSMark Adams 
1480c1c463dbSMark Adams   /* general events */
1481fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1482fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1483fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1484fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1485c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1486c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1487fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1488c1c463dbSMark Adams 
14895b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
14905b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
14915b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
14925b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
14935b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
14945b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
14955b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
14965b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
14975b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1498bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
14995b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
15005b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
15015b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
15025b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
15035b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
15045b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
15055b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
15065b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
15075b89ad90SMark F. Adams 
15085b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
15095b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
15105b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
15115b89ad90SMark F. Adams   /* create timer stages */
15125b89ad90SMark F. Adams #if defined GAMG_STAGES
15135b89ad90SMark F. Adams   {
15145b89ad90SMark F. Adams     char     str[32];
15155b89ad90SMark F. Adams     PetscInt lidx;
15165b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
15175b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
15185b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
15195b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
15205b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
15215b89ad90SMark F. Adams     }
15225b89ad90SMark F. Adams   }
15235b89ad90SMark F. Adams #endif
15245b89ad90SMark F. Adams #endif
15253e3471ccSMark Adams   PetscFunctionReturn(0);
15263e3471ccSMark Adams }
15273e3471ccSMark Adams 
15283e3471ccSMark Adams /*@C
15291c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
15301c1aac46SBarry Smith     called from PetscFinalize() automatically.
15313e3471ccSMark Adams 
15323e3471ccSMark Adams  Level: developer
15333e3471ccSMark Adams 
15343e3471ccSMark Adams  .seealso: PetscFinalize()
15353e3471ccSMark Adams @*/
15363e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
15373e3471ccSMark Adams {
15383e3471ccSMark Adams   PetscErrorCode ierr;
15393e3471ccSMark Adams 
15403e3471ccSMark Adams   PetscFunctionBegin;
15413e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
15423e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
15433e3471ccSMark Adams   PetscFunctionReturn(0);
15443e3471ccSMark Adams }
1545a36cf38bSToby Isaac 
1546a36cf38bSToby Isaac /*@C
1547a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1548a36cf38bSToby Isaac 
1549a36cf38bSToby Isaac  Input Parameters:
1550a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1551a36cf38bSToby Isaac  - create - function for creating the gamg context.
1552a36cf38bSToby Isaac 
1553a36cf38bSToby Isaac   Level: advanced
1554a36cf38bSToby Isaac 
15551c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1556a36cf38bSToby Isaac @*/
1557a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1558a36cf38bSToby Isaac {
1559a36cf38bSToby Isaac   PetscErrorCode ierr;
1560a36cf38bSToby Isaac 
1561a36cf38bSToby Isaac   PetscFunctionBegin;
1562a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1563a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1564a36cf38bSToby Isaac   PetscFunctionReturn(0);
1565a36cf38bSToby Isaac }
1566a36cf38bSToby Isaac 
1567