xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision f7df55f016625f8a883057482344e9e0d655e2a7)
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 */
718c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>    /*I "petscksp.h" I*/
8f96513f1SMatthew G Knepley 
90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
11b4fbaa2aSMark F. Adams #endif
120cbbd2e1SMark F. Adams 
130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
15fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
20fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
210cbbd2e1SMark F. Adams #endif
220cbbd2e1SMark F. Adams 
23b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
240cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
2518c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
26b4fbaa2aSMark F. Adams #endif
27f96513f1SMatthew G Knepley 
280a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
293e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
309d5b6da9SMark F. Adams 
31d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
32d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
33d3d6bff4SMark F. Adams {
3418c3aa7eSMark   PetscErrorCode ierr, level;
35d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
36d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
37d3d6bff4SMark F. Adams 
38d3d6bff4SMark F. Adams   PetscFunctionBegin;
3922a233eaSStefano Zampini   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
401c1aac46SBarry Smith   pc_gamg->data_sz = 0;
41878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
4218c3aa7eSMark   for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
4318c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
4418c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
4518c3aa7eSMark   }
4618c3aa7eSMark   pc_gamg->emin = 0;
4718c3aa7eSMark   pc_gamg->emax = 0;
48a2f3521dSMark F. Adams   PetscFunctionReturn(0);
49a2f3521dSMark F. Adams }
50a2f3521dSMark F. Adams 
515b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
525b89ad90SMark F. Adams /*
53c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
54a147abb0SMark F. Adams      of active processors.
555b89ad90SMark F. Adams 
565b89ad90SMark F. Adams    Input Parameter:
57a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
58a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
599d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
60c5bfad50SMark F. Adams    . cr_bs - coarse block size
613530afc2SMark F. Adams    In/Output Parameter:
62a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
63afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6411e60469SMark F. Adams    Output Parameter:
653530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
665b89ad90SMark F. Adams */
675cb416c2SMark F. Adams 
68171cca9aSMark 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)
695b89ad90SMark F. Adams {
70a2f3521dSMark F. Adams   PetscErrorCode  ierr;
719d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
72486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
73a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
743b4367a7SBarry Smith   MPI_Comm        comm;
75c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
763ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
775b89ad90SMark F. Adams 
785b89ad90SMark F. Adams   PetscFunctionBegin;
793b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
803b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
813b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
82c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
839d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
84038e3b61SMark F. Adams 
85ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
86ce7c7f2fSMark Adams 
873ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
880298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
893ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
903ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9173911c69SBarry Smith   } else {
923ae0bb68SMark Adams     PetscInt  bs;
933ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
943ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
953ae0bb68SMark Adams   }
96c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
97171cca9aSMark Adams   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
98171cca9aSMark Adams   else {
99472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
1000298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
101a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1027f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
103c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
104a2f3521dSMark F. Adams   }
105f852f58cSMark F. Adams 
1062e3501ffSMark Adams   if (new_size==nactive) {
107ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
108ce7c7f2fSMark Adams     if (new_size < size) {
109ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
11031cb4603SMark Adams       ierr = PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr);
111ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
112b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
113b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
114ce7c7f2fSMark Adams       }
115ce7c7f2fSMark Adams     }
116ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1172e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
118192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
119885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
12071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
12171959b99SBarry 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);
122ce7c7f2fSMark Adams     /* get new_size and rfactor */
123ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
124ce7c7f2fSMark Adams       /* find factor */
125ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
126ce7c7f2fSMark Adams       else {
127ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
128ce7c7f2fSMark Adams         jj = -1;
129ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
130ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
131ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
132ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
133ce7c7f2fSMark Adams             if (fact > best_fact) {
134ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
135ce7c7f2fSMark Adams             }
136ce7c7f2fSMark Adams           }
137ce7c7f2fSMark Adams         }
138ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
139ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
140ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
141ce7c7f2fSMark Adams         else expand_factor = rfactor;
142ce7c7f2fSMark Adams       }
143ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
1444cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1454cdfd227SMark         *a_Amat_crs = Cmat;
14631cb4603SMark Adams         ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
147ce7c7f2fSMark Adams         PetscFunctionReturn(0);
148ce7c7f2fSMark Adams       }
149ce7c7f2fSMark Adams     }
1504cdfd227SMark #if defined PETSC_GAMG_USE_LOG
1514cdfd227SMark     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
1524cdfd227SMark #endif
153a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
154785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1552e3501ffSMark Adams     if (pc_gamg->repart) {
156a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
1575a9b9e01SMark F. Adams       Mat      adj;
15831cb4603SMark 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);
159a2f3521dSMark F. Adams       /* get 'adj' */
160c5bfad50SMark F. Adams       if (cr_bs == 1) {
161038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
162806fa848SBarry Smith       } else {
163a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
164eb07cef2SMark F. Adams         Mat               tMat;
165a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
166b4fbaa2aSMark F. Adams         const PetscScalar *vals;
167b4fbaa2aSMark F. Adams         const PetscInt    *idx;
168a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
16939d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
170d9558ea9SBarry Smith         MatType           mtype;
171b4fbaa2aSMark F. Adams 
172e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
173a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
174a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
175c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
1760a545947SLisandro Dalcin           ierr      = MatGetRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
177c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
178c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
1790a545947SLisandro Dalcin           ierr      = MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
1803ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
1813ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
18258471d46SMark F. Adams         }
1836876a03eSMark F. Adams 
184d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
1853b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
1863ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
187d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
188a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
189a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
190e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
191eb07cef2SMark F. Adams 
192a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
193c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
19422063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
195eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
196c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
197eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
198eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
199eb07cef2SMark F. Adams           }
20022063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
201eb07cef2SMark F. Adams         }
202eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
203eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
204eb07cef2SMark F. Adams 
205b4fbaa2aSMark F. Adams         if (llev++ == -1) {
206b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2078caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2083b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
209b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2103bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
211b4fbaa2aSMark F. Adams         }
212eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
213eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
214a2f3521dSMark F. Adams       } /* create 'adj' */
215f150b916SMark F. Adams 
216a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2175a9b9e01SMark F. Adams         char            prefix[256];
2185a9b9e01SMark F. Adams         const char      *pcpre;
219b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
220b4fbaa2aSMark F. Adams         MatPartitioning mpart;
221a4b7d37bSMark F. Adams         IS              proc_is;
2222f03bc48SMark F. Adams 
2233b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2245ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2259d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2268caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
22759a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
22811e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
229c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
230a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
23111e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2325a9b9e01SMark F. Adams 
2335ef31b24SMark F. Adams         /* collect IS info */
234785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
235a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
236a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
237c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
238ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
239eb07cef2SMark F. Adams           }
2405ef31b24SMark F. Adams         }
241a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
242a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2435ef31b24SMark F. Adams       }
2445ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2455a9b9e01SMark F. Adams 
2463b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2478263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
24831cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
249ce7c7f2fSMark Adams       PetscInt targetPE;
2504cdfd227SMark       if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
251302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
252ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2533b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
254a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
255e33ef3b1SMark F. Adams 
25611e60469SMark F. Adams     /*
257a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
25811e60469SMark F. Adams     */
259a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
2607700e67bSMark Adams     is_eq_num_prim = is_eq_num;
26111e60469SMark F. Adams     /*
262a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
26311e60469SMark F. Adams     */
264c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
265c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
266a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
267ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
268a2f3521dSMark F. Adams 
269a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
270885364a3SMark 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 */
271885364a3SMark Adams     {
272885364a3SMark Adams       Vec            src_crd, dest_crd;
273885364a3SMark 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;
274885364a3SMark Adams       VecScatter     vecscat;
275885364a3SMark Adams       PetscScalar    *array;
276885364a3SMark Adams       IS isscat;
277a2f3521dSMark F. Adams       /* move data (for primal equations only) */
27822063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
2793b4367a7SBarry Smith       ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
2803ae0bb68SMark Adams       ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
281c0dedaeaSBarry Smith       ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
28211e60469SMark F. Adams       /*
2839d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
284c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
28511e60469SMark F. Adams       */
286854ce69bSBarry Smith       ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
287a2f3521dSMark F. Adams       ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2883ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
289c5bfad50SMark F. Adams         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
290a2f3521dSMark F. Adams         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
29111e60469SMark F. Adams       }
292a2f3521dSMark F. Adams       ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
2933ae0bb68SMark Adams       ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
29492a756f0SMark F. Adams       ierr = PetscFree(tidx);CHKERRQ(ierr);
29511e60469SMark F. Adams       /*
29611e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
29711e60469SMark F. Adams       */
2983ae0bb68SMark Adams       ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
2999d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3003ae0bb68SMark Adams         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3013ae0bb68SMark Adams         for (ii=0; ii<ncrs; ii++) {
3029d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
303a2f3521dSMark F. Adams             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
304c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
305676e1743SMark F. Adams             ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
306d3d6bff4SMark F. Adams           }
307038e3b61SMark F. Adams         }
308eb07cef2SMark F. Adams       }
309eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
310eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
31111e60469SMark F. Adams       /*
31211e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
31311e60469SMark F. Adams         to the correct processor
31411e60469SMark F. Adams       */
3159448b7f1SJunchao Zhang       ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
31611e60469SMark F. Adams       ierr = ISDestroy(&isscat);CHKERRQ(ierr);
31711e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31811e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31911e60469SMark F. Adams       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
32011e60469SMark F. Adams       ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
32111e60469SMark F. Adams       /*
32211e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
32311e60469SMark F. Adams       */
324c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
325578f55a3SPeter Brune       ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3262fa5cd67SKarl Rupp 
3273ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3283ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3292fa5cd67SKarl Rupp 
33011e60469SMark F. Adams       ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3319d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3323ae0bb68SMark Adams         for (ii=0; ii<ncrs_new; ii++) {
3339d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
334a2f3521dSMark F. Adams             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
335c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
336d3d6bff4SMark F. Adams           }
337038e3b61SMark F. Adams         }
338038e3b61SMark F. Adams       }
33911e60469SMark F. Adams       ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
34011e60469SMark F. Adams       ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
341885364a3SMark Adams     }
342a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3440cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
345ed3f9983SMark F. Adams #endif
34611e60469SMark F. Adams     /*
3477dae84e0SHong Zhang       Invert for MatCreateSubMatrix
34811e60469SMark F. Adams     */
349a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
350a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
351c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
352a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
353a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
354a2f3521dSMark F. Adams     }
3553cb8563fSToby Isaac     if (Pcolumnperm) {
3563cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3573cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3583cb8563fSToby Isaac     }
359a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3610cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3620cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
363ed3f9983SMark F. Adams #endif
364a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
365a2f3521dSMark F. Adams     {
366a2f3521dSMark F. Adams       Mat mat;
3677dae84e0SHong Zhang       ierr        = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
368a2f3521dSMark F. Adams       *a_Amat_crs = mat;
369a2f3521dSMark F. Adams     }
370038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
371a2f3521dSMark F. Adams 
3720cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3730cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
374ed3f9983SMark F. Adams #endif
37511e60469SMark F. Adams     /* prolongator */
37611e60469SMark F. Adams     {
37711e60469SMark F. Adams       IS       findices;
378a2f3521dSMark F. Adams       PetscInt Istart,Iend;
379a2f3521dSMark F. Adams       Mat      Pnew;
38062294041SBarry Smith 
381a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
3820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3830cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
384ed3f9983SMark F. Adams #endif
3853b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
386c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
3877dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
38811e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
389c5bfad50SMark F. Adams 
3900cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3910cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
392ed3f9983SMark F. Adams #endif
3933530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
394a2f3521dSMark F. Adams 
395a2f3521dSMark F. Adams       /* output - repartitioned */
396a2f3521dSMark F. Adams       *a_P_inout = Pnew;
397e33ef3b1SMark F. Adams     }
398a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
3995b89ad90SMark F. Adams 
400c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
401ce7c7f2fSMark Adams 
402ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
403ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
404ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4058bca76a6SMark Adams       static PetscInt llev = 2;
40639d09545SMark Adams       ierr = PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr);
407ce7c7f2fSMark Adams #endif
408b470e4b4SRichard Tran Mills       ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
409b470e4b4SRichard Tran Mills       ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
410ce7c7f2fSMark Adams       if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
411ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
412ce7c7f2fSMark Adams         PetscMPIInt size;
413ce7c7f2fSMark Adams         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
414ce7c7f2fSMark Adams         if (size > 1) {
415ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
416b470e4b4SRichard Tran Mills           ierr = VecBindToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr);
417b470e4b4SRichard Tran Mills           ierr = VecBindToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr);
418ce7c7f2fSMark Adams         }
419ce7c7f2fSMark Adams       }
420ce7c7f2fSMark Adams     }
4214cdfd227SMark #if defined PETSC_GAMG_USE_LOG
4224cdfd227SMark     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
4234cdfd227SMark #endif
424a2f3521dSMark F. Adams   }
4255b89ad90SMark F. Adams   PetscFunctionReturn(0);
4265b89ad90SMark F. Adams }
4275b89ad90SMark F. Adams 
4285b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4295b89ad90SMark F. Adams /*
4305b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4315b89ad90SMark F. Adams                     by setting data structures and options.
4325b89ad90SMark F. Adams 
4335b89ad90SMark F. Adams    Input Parameter:
4345b89ad90SMark F. Adams .  pc - the preconditioner context
4355b89ad90SMark F. Adams 
4365b89ad90SMark F. Adams */
4379d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4385b89ad90SMark F. Adams {
4395b89ad90SMark F. Adams   PetscErrorCode ierr;
4409d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4415b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4422adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
44318c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
4443b4367a7SBarry Smith   MPI_Comm       comm;
445c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
44618c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
44718c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
448a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
449569f4572SMark Adams   MatInfo        info;
450171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4515ef31b24SMark F. Adams 
4525b89ad90SMark F. Adams   PetscFunctionBegin;
4533b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4543b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4553b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
456dfd5c07aSMark F. Adams 
45784d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4581c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
459878e152fSMark F. Adams       /* reset everything */
460878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
461878e152fSMark F. Adams       pc->setupcalled = 0;
462806fa848SBarry Smith     } else {
46384d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
46403a628feSMark F. Adams       /* just do Galerkin grids */
46558471d46SMark F. Adams       Mat          B,dA,dB;
46658471d46SMark F. Adams 
46771959b99SBarry Smith       if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4689d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
46958471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
47023ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
47158471d46SMark F. Adams         /* (re)set to get dirty flag */
47223ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
47358471d46SMark F. Adams 
4742fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
475ef3f0257SMark Adams           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
476ef3f0257SMark Adams           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
477ef3f0257SMark Adams             ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
478ef3f0257SMark Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr);
479084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
48003a628feSMark F. Adams             mglevels[level]->A = B;
481806fa848SBarry Smith           } else {
482ef3f0257SMark Adams             ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
48323ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
48458471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
48503a628feSMark F. Adams           }
48623ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
48758471d46SMark F. Adams           dB   = B;
48858471d46SMark F. Adams         }
4895f8cf99dSMark F. Adams       }
490d5280255SMark F. Adams 
491d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
49258471d46SMark F. Adams       PetscFunctionReturn(0);
493eb07cef2SMark F. Adams     }
494878e152fSMark F. Adams   }
495f6536408SMark F. Adams 
496878e152fSMark F. Adams   if (!pc_gamg->data) {
497878e152fSMark F. Adams     if (pc_gamg->orig_data) {
498878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
4990298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5002fa5cd67SKarl Rupp 
501878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
502878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
503878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5042fa5cd67SKarl Rupp 
505785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
506878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
507806fa848SBarry Smith     } else {
5081ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5097700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5109d5b6da9SMark F. Adams     }
511878e152fSMark F. Adams   }
512878e152fSMark F. Adams 
513878e152fSMark F. Adams   /* cache original data for reuse */
5141c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
515785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
516878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
517878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
518878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
519878e152fSMark F. Adams   }
520038e3b61SMark F. Adams 
521302f38e8SMark F. Adams   /* get basic dims */
522302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
523171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
52484d3f75bSMark F. Adams 
525569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
526569f4572SMark Adams   nnz0   = info.nz_used;
527569f4572SMark Adams   nnztot = info.nz_used;
52862294041SBarry 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);
529569f4572SMark Adams 
530a2f3521dSMark F. Adams   /* Get A_i and R_i */
53162294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5329ab59c8bSMark Adams     pc_gamg->current_level = level;
53318c3aa7eSMark     if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
5345b89ad90SMark F. Adams     level1 = level + 1;
5350cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5360cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
537a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
538a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
539b4fbaa2aSMark F. Adams #endif
540a2f3521dSMark F. Adams #endif
541c8b0795cSMark F. Adams     { /* construct prolongator */
542725b86d8SJed Brown       Mat              Gmat;
5430cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5447700e67bSMark Adams       Mat              Prol11;
545c8b0795cSMark F. Adams 
5467700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5471ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5487700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
549c8b0795cSMark F. Adams 
550a2f3521dSMark F. Adams       /* could have failed to create new level */
551a2f3521dSMark F. Adams       if (Prol11) {
552*f7df55f0SStefano Zampini         const char *prefix;
553*f7df55f0SStefano Zampini         char       addp[32];
554*f7df55f0SStefano Zampini 
5559d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5560298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
557a2f3521dSMark F. Adams 
558fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
559c8b0795cSMark F. Adams           /* smooth */
560fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
561c8b0795cSMark F. Adams         }
562c8b0795cSMark F. Adams 
5630c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
5641b18a24aSMark Adams           PetscInt bs;
5651b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5660a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
567ffc955d6SMark F. Adams         }
568ffc955d6SMark F. Adams 
569*f7df55f0SStefano Zampini         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
570*f7df55f0SStefano Zampini         ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr);
571*f7df55f0SStefano Zampini         ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",level);CHKERRQ(ierr);
572*f7df55f0SStefano Zampini         ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr);
573*f7df55f0SStefano Zampini         /* if we reuse the prolongator, then it is better to generate the transpose by default with CUDA
574*f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
575*f7df55f0SStefano Zampini #if defined(PETSC_HAVE_CUDA)
576*f7df55f0SStefano Zampini         {
577*f7df55f0SStefano Zampini           PetscBool ismpiaij;
578*f7df55f0SStefano Zampini 
579*f7df55f0SStefano Zampini           ierr = PetscObjectBaseTypeCompare((PetscObject)Prol11,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
580*f7df55f0SStefano Zampini           if (ismpiaij) {
581*f7df55f0SStefano Zampini             Mat Prol_d,Prol_o;
582*f7df55f0SStefano Zampini 
583*f7df55f0SStefano Zampini             ierr = MatMPIAIJGetSeqAIJ(Prol11,&Prol_d,&Prol_o,NULL);CHKERRQ(ierr);
584*f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol_d,pc_gamg->reuse_prol);CHKERRQ(ierr);
585*f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol_o,pc_gamg->reuse_prol);CHKERRQ(ierr);
586*f7df55f0SStefano Zampini           } else {
587*f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol11,pc_gamg->reuse_prol);CHKERRQ(ierr);
588*f7df55f0SStefano Zampini           }
589*f7df55f0SStefano Zampini         }
590*f7df55f0SStefano Zampini #endif
591*f7df55f0SStefano Zampini         ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr);
5924bde40a0SMark Adams         Parr[level1] = Prol11;
5934bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
5944bde40a0SMark Adams 
595a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
59641b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
597a2f3521dSMark F. Adams     } /* construct prolongator scope */
5980cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5990cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
600c8b0795cSMark F. Adams #endif
6017f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
602171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
603569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
60462294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
605a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
606a90e85d9SMark Adams #endif
607c8b0795cSMark F. Adams       break;
608c8b0795cSMark F. Adams     }
6090cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6100cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
611b4fbaa2aSMark F. Adams #endif
612171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
613171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
614171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6150e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
616171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
617a2f3521dSMark F. Adams 
6180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6190cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
620b4fbaa2aSMark F. Adams #endif
621171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
622569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
623569f4572SMark Adams     nnztot += info.nz_used;
6241d5b2942SMark 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);
625569f4572SMark Adams 
6260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
627b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
628b4fbaa2aSMark F. Adams #endif
629a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6309ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
6319ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
632a90e85d9SMark Adams       level++;
633a90e85d9SMark Adams       break;
634a90e85d9SMark Adams     }
635c8b0795cSMark F. Adams   } /* levels */
636c8b0795cSMark F. Adams   ierr                  = PetscFree(pc_gamg->data);CHKERRQ(ierr);
637c8b0795cSMark F. Adams 
638569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
6399d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6405b89ad90SMark F. Adams   fine_level       = level;
6410298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6425b89ad90SMark F. Adams 
64362294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
644d5280255SMark F. Adams     /* set default smoothers & set operators */
64562294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
646ffc955d6SMark F. Adams       KSP smoother;
647ffc955d6SMark F. Adams       PC  subpc;
648a2f3521dSMark F. Adams 
6499d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
650f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
651ffc955d6SMark F. Adams 
652a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
653a2f3521dSMark F. Adams       /* set ops */
65423ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
655a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
656a2f3521dSMark F. Adams 
657a2f3521dSMark F. Adams       /* set defaults */
6586c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
659a2f3521dSMark F. Adams 
6600c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6610c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6622d3561bbSSatish Balay         PetscInt sz;
6637a28f3e5SMark Adams         IS       *iss;
664a2f3521dSMark F. Adams 
6652d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6667a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
6670c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6680a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6690c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6707f66b68fSMark Adams         if (!sz) {
671ffc955d6SMark F. Adams           IS       is;
6720a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6737a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
674a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
675806fa848SBarry Smith         } else {
676a94c3b12SMark F. Adams           PetscInt kk;
6777a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
678a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
6797a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
680a94c3b12SMark F. Adams           }
6817a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
682ffc955d6SMark F. Adams         }
6830298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
684ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
685806fa848SBarry Smith       } else {
686890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
687ffc955d6SMark F. Adams       }
688d5280255SMark F. Adams     }
689d5280255SMark F. Adams     {
690d5280255SMark F. Adams       /* coarse grid */
691d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
692d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
693d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
69423ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
695cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
696d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
697d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
698d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
699d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
70071959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
70171959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
702d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
703d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7049dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7052fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
70608e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
7075b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7085b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7095b42dca8SJed Brown          * view every subdomain as though they were different. */
7105b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
711d5280255SMark F. Adams       }
712cf8ae1d3SMark Adams     }
713d5280255SMark F. Adams 
714d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
715d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
716e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
717d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
71869aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
719d5280255SMark F. Adams 
72018c3aa7eSMark     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
72118c3aa7eSMark 
72218c3aa7eSMark     /* setup cheby eigen estimates from SA */
72318c3aa7eSMark     for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
72418c3aa7eSMark       KSP       smoother;
72518c3aa7eSMark       PetscBool ischeb;
72618c3aa7eSMark       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
72718c3aa7eSMark       ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr);
72818c3aa7eSMark       if (ischeb) {
72918c3aa7eSMark         KSP_Chebyshev  *cheb = (KSP_Chebyshev*)smoother->data;
73018c3aa7eSMark         if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) { /* let command line emax override using SA's eigenvalues */
73118c3aa7eSMark           PC        subpc;
73218c3aa7eSMark           PetscBool isjac;
73318c3aa7eSMark           ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
73418c3aa7eSMark           ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr);
73518c3aa7eSMark           if ((isjac && pc_gamg->use_sa_esteig==-1) || pc_gamg->use_sa_esteig==1) {
73618c3aa7eSMark             PetscReal       emax,emin;
73718c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
73818c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
73918c3aa7eSMark             ierr = PetscInfo4(pc,"PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",level,Aarr[level]->rmap->N,(double)emax,(double)emin);CHKERRQ(ierr);
74018c3aa7eSMark             cheb->emin_computed = emin;
74118c3aa7eSMark             cheb->emax_computed = emax;
74218c3aa7eSMark             ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr);
74318c3aa7eSMark           }
74418c3aa7eSMark         }
74518c3aa7eSMark       }
74618c3aa7eSMark     }
74718c3aa7eSMark 
748d5280255SMark F. Adams     /* clean up */
749d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
750587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
751587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
7525b89ad90SMark F. Adams     }
75318c3aa7eSMark 
754806fa848SBarry Smith   } else {
7555f8cf99dSMark F. Adams     KSP smoother;
756302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
7579d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
75823ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
7595f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
7609d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
7615f8cf99dSMark F. Adams   }
7625b89ad90SMark F. Adams   PetscFunctionReturn(0);
7635b89ad90SMark F. Adams }
7645b89ad90SMark F. Adams 
765eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
7665b89ad90SMark F. Adams /*
7675b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
7685b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
7695b89ad90SMark F. Adams 
7705b89ad90SMark F. Adams    Input Parameter:
7715b89ad90SMark F. Adams .  pc - the preconditioner context
7725b89ad90SMark F. Adams 
7735b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
7745b89ad90SMark F. Adams */
7755b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
7765b89ad90SMark F. Adams {
7775b89ad90SMark F. Adams   PetscErrorCode ierr;
7785b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
7795b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
7805b89ad90SMark F. Adams 
7815b89ad90SMark F. Adams   PetscFunctionBegin;
7825b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
7839b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
7849b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
7859b8ffb57SJed Brown   }
7861ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
7871ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
7885b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
7895b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
7905b89ad90SMark F. Adams   PetscFunctionReturn(0);
7915b89ad90SMark F. Adams }
7925b89ad90SMark F. Adams 
793676e1743SMark F. Adams /*@
794cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
795676e1743SMark F. Adams 
7961cc46a46SBarry Smith    Logically Collective on PC
797676e1743SMark F. Adams 
798676e1743SMark F. Adams    Input Parameters:
7991cc46a46SBarry Smith +  pc - the preconditioner context
8001cc46a46SBarry Smith -  n - the number of equations
801676e1743SMark F. Adams 
802676e1743SMark F. Adams 
803676e1743SMark F. Adams    Options Database Key:
8041cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
805676e1743SMark F. Adams 
80695452b02SPatrick Sanan    Notes:
80795452b02SPatrick 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
808cab9ed1eSBarry Smith           that has degrees of freedom
809cab9ed1eSBarry Smith 
810676e1743SMark F. Adams    Level: intermediate
811676e1743SMark F. Adams 
8121c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
813676e1743SMark F. Adams @*/
814676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
815676e1743SMark F. Adams {
816676e1743SMark F. Adams   PetscErrorCode ierr;
817676e1743SMark F. Adams 
818676e1743SMark F. Adams   PetscFunctionBegin;
819676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
820676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
821676e1743SMark F. Adams   PetscFunctionReturn(0);
822676e1743SMark F. Adams }
823676e1743SMark F. Adams 
8241e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
825676e1743SMark F. Adams {
826c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
827c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
828676e1743SMark F. Adams 
829676e1743SMark F. Adams   PetscFunctionBegin;
8309d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
831676e1743SMark F. Adams   PetscFunctionReturn(0);
832676e1743SMark F. Adams }
833676e1743SMark F. Adams 
834389730f3SMark F. Adams /*@
835cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
836389730f3SMark F. Adams 
837389730f3SMark F. Adams  Collective on PC
838389730f3SMark F. Adams 
839389730f3SMark F. Adams    Input Parameters:
8401cc46a46SBarry Smith +  pc - the preconditioner context
8411cc46a46SBarry Smith -  n - maximum number of equations to aim for
842389730f3SMark F. Adams 
843389730f3SMark F. Adams    Options Database Key:
8441cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
845389730f3SMark F. Adams 
84674329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
84774329af1SBarry Smith      has less than 1000 unknowns.
84874329af1SBarry Smith 
849389730f3SMark F. Adams    Level: intermediate
850389730f3SMark F. Adams 
8511c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
852389730f3SMark F. Adams @*/
853389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
854389730f3SMark F. Adams {
855389730f3SMark F. Adams   PetscErrorCode ierr;
856389730f3SMark F. Adams 
857389730f3SMark F. Adams   PetscFunctionBegin;
858389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
859389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
860389730f3SMark F. Adams   PetscFunctionReturn(0);
861389730f3SMark F. Adams }
862389730f3SMark F. Adams 
8631e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
864389730f3SMark F. Adams {
865389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
866389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
867389730f3SMark F. Adams 
868389730f3SMark F. Adams   PetscFunctionBegin;
8699d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
870389730f3SMark F. Adams   PetscFunctionReturn(0);
871389730f3SMark F. Adams }
872389730f3SMark F. Adams 
873676e1743SMark F. Adams /*@
874cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
875676e1743SMark F. Adams 
876676e1743SMark F. Adams    Collective on PC
877676e1743SMark F. Adams 
878676e1743SMark F. Adams    Input Parameters:
8791cc46a46SBarry Smith +  pc - the preconditioner context
8801cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
881676e1743SMark F. Adams 
882676e1743SMark F. Adams    Options Database Key:
8831cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
884676e1743SMark F. Adams 
88595452b02SPatrick Sanan    Notes:
88695452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
887cab9ed1eSBarry Smith 
888676e1743SMark F. Adams    Level: intermediate
889676e1743SMark F. Adams 
890676e1743SMark F. Adams .seealso: ()
891676e1743SMark F. Adams @*/
892cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
893676e1743SMark F. Adams {
894676e1743SMark F. Adams   PetscErrorCode ierr;
895676e1743SMark F. Adams 
896676e1743SMark F. Adams   PetscFunctionBegin;
897676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
898cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
899676e1743SMark F. Adams   PetscFunctionReturn(0);
900676e1743SMark F. Adams }
901676e1743SMark F. Adams 
902cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
903676e1743SMark F. Adams {
904c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
905c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
906676e1743SMark F. Adams 
907676e1743SMark F. Adams   PetscFunctionBegin;
9089d5b6da9SMark F. Adams   pc_gamg->repart = n;
909676e1743SMark F. Adams   PetscFunctionReturn(0);
910676e1743SMark F. Adams }
911676e1743SMark F. Adams 
912dfd5c07aSMark F. Adams /*@
91318c3aa7eSMark    PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby)
91418c3aa7eSMark 
91518c3aa7eSMark    Collective on PC
91618c3aa7eSMark 
91718c3aa7eSMark    Input Parameters:
91818c3aa7eSMark +  pc - the preconditioner context
91918c3aa7eSMark -  n - number of its
92018c3aa7eSMark 
92118c3aa7eSMark    Options Database Key:
92218c3aa7eSMark .  -pc_gamg_esteig_ksp_max_it <its>
92318c3aa7eSMark 
92418c3aa7eSMark    Notes:
92518c3aa7eSMark 
92618c3aa7eSMark    Level: intermediate
92718c3aa7eSMark 
92818c3aa7eSMark .seealso: ()
92918c3aa7eSMark @*/
93018c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n)
93118c3aa7eSMark {
93218c3aa7eSMark   PetscErrorCode ierr;
93318c3aa7eSMark 
93418c3aa7eSMark   PetscFunctionBegin;
93518c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
93618c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
93718c3aa7eSMark   PetscFunctionReturn(0);
93818c3aa7eSMark }
93918c3aa7eSMark 
94018c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n)
94118c3aa7eSMark {
94218c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
94318c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
94418c3aa7eSMark 
94518c3aa7eSMark   PetscFunctionBegin;
94618c3aa7eSMark   pc_gamg->esteig_max_it = n;
94718c3aa7eSMark   PetscFunctionReturn(0);
94818c3aa7eSMark }
94918c3aa7eSMark 
95018c3aa7eSMark /*@
95118c3aa7eSMark    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother
95218c3aa7eSMark 
95318c3aa7eSMark    Collective on PC
95418c3aa7eSMark 
95518c3aa7eSMark    Input Parameters:
95618c3aa7eSMark +  pc - the preconditioner context
95718c3aa7eSMark -  n - number of its
95818c3aa7eSMark 
95918c3aa7eSMark    Options Database Key:
96018c3aa7eSMark .  -pc_gamg_use_sa_esteig <true,false>
96118c3aa7eSMark 
96218c3aa7eSMark    Notes:
96318c3aa7eSMark 
96418c3aa7eSMark    Level: intermediate
96518c3aa7eSMark 
96618c3aa7eSMark .seealso: ()
96718c3aa7eSMark @*/
96818c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
96918c3aa7eSMark {
97018c3aa7eSMark   PetscErrorCode ierr;
97118c3aa7eSMark 
97218c3aa7eSMark   PetscFunctionBegin;
97318c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
97418c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
97518c3aa7eSMark   PetscFunctionReturn(0);
97618c3aa7eSMark }
97718c3aa7eSMark 
97818c3aa7eSMark static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscInt n)
97918c3aa7eSMark {
98018c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
98118c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
98218c3aa7eSMark 
98318c3aa7eSMark   PetscFunctionBegin;
98418c3aa7eSMark   pc_gamg->use_sa_esteig = n ? 1 : 0;
98518c3aa7eSMark   PetscFunctionReturn(0);
98618c3aa7eSMark }
98718c3aa7eSMark 
98818c3aa7eSMark /*@C
98918c3aa7eSMark    PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby)
99018c3aa7eSMark 
99118c3aa7eSMark    Collective on PC
99218c3aa7eSMark 
99318c3aa7eSMark    Input Parameters:
99418c3aa7eSMark +  pc - the preconditioner context
99518c3aa7eSMark -  t - ksp type
99618c3aa7eSMark 
99718c3aa7eSMark    Options Database Key:
99818c3aa7eSMark .  -pc_gamg_esteig_ksp_type <type>
99918c3aa7eSMark 
100018c3aa7eSMark    Notes:
100118c3aa7eSMark 
100218c3aa7eSMark    Level: intermediate
100318c3aa7eSMark 
100418c3aa7eSMark .seealso: ()
100518c3aa7eSMark @*/
100618c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[])
100718c3aa7eSMark {
100818c3aa7eSMark   PetscErrorCode ierr;
100918c3aa7eSMark 
101018c3aa7eSMark   PetscFunctionBegin;
101118c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
101218c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr);
101318c3aa7eSMark   PetscFunctionReturn(0);
101418c3aa7eSMark }
101518c3aa7eSMark 
101618c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[])
101718c3aa7eSMark {
101818c3aa7eSMark   PetscErrorCode ierr;
101918c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
102018c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
102118c3aa7eSMark 
102218c3aa7eSMark   PetscFunctionBegin;
102318c3aa7eSMark   ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr);
102418c3aa7eSMark   PetscFunctionReturn(0);
102518c3aa7eSMark }
102618c3aa7eSMark 
102718c3aa7eSMark /*@
102818c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
102918c3aa7eSMark 
103018c3aa7eSMark    Collective on PC
103118c3aa7eSMark 
103218c3aa7eSMark    Input Parameters:
103318c3aa7eSMark +  pc - the preconditioner context
103418c3aa7eSMark -  emax - max eigenvalue
103518c3aa7eSMark .  emin - min eigenvalue
103618c3aa7eSMark 
103718c3aa7eSMark    Options Database Key:
103818c3aa7eSMark .  -pc_gamg_eigenvalues
103918c3aa7eSMark 
104018c3aa7eSMark    Level: intermediate
104118c3aa7eSMark 
104218c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType()
104318c3aa7eSMark @*/
104418c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
104518c3aa7eSMark {
104618c3aa7eSMark   PetscErrorCode ierr;
104718c3aa7eSMark 
104818c3aa7eSMark   PetscFunctionBegin;
104918c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
105018c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr);
105118c3aa7eSMark   PetscFunctionReturn(0);
105218c3aa7eSMark }
105318c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
105418c3aa7eSMark {
105518c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
105618c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
105718c3aa7eSMark 
105818c3aa7eSMark   PetscFunctionBegin;
105918c3aa7eSMark   if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
106018c3aa7eSMark   if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
106118c3aa7eSMark   pc_gamg->emax = emax;
106218c3aa7eSMark   pc_gamg->emin = emin;
106318c3aa7eSMark 
106418c3aa7eSMark   PetscFunctionReturn(0);
106518c3aa7eSMark }
106618c3aa7eSMark 
106718c3aa7eSMark /*@
1068cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1069dfd5c07aSMark F. Adams 
1070dfd5c07aSMark F. Adams    Collective on PC
1071dfd5c07aSMark F. Adams 
1072dfd5c07aSMark F. Adams    Input Parameters:
10731cc46a46SBarry Smith +  pc - the preconditioner context
10741cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1075dfd5c07aSMark F. Adams 
1076dfd5c07aSMark F. Adams    Options Database Key:
10771cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1078dfd5c07aSMark F. Adams 
1079dfd5c07aSMark F. Adams    Level: intermediate
1080dfd5c07aSMark F. Adams 
108195452b02SPatrick Sanan    Notes:
108295452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1083cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
1084cab9ed1eSBarry Smith 
1085dfd5c07aSMark F. Adams .seealso: ()
1086dfd5c07aSMark F. Adams @*/
10871cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1088dfd5c07aSMark F. Adams {
1089dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1090dfd5c07aSMark F. Adams 
1091dfd5c07aSMark F. Adams   PetscFunctionBegin;
1092dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10931cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1094dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1095dfd5c07aSMark F. Adams }
1096dfd5c07aSMark F. Adams 
10971cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1098dfd5c07aSMark F. Adams {
1099dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1100dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1101dfd5c07aSMark F. Adams 
1102dfd5c07aSMark F. Adams   PetscFunctionBegin;
1103dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1104dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1105dfd5c07aSMark F. Adams }
1106dfd5c07aSMark F. Adams 
1107ffc955d6SMark F. Adams /*@
1108cab9ed1eSBarry 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.
1109ffc955d6SMark F. Adams 
1110ffc955d6SMark F. Adams    Collective on PC
1111ffc955d6SMark F. Adams 
1112ffc955d6SMark F. Adams    Input Parameters:
1113cab9ed1eSBarry Smith +  pc - the preconditioner context
1114cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1115ffc955d6SMark F. Adams 
1116ffc955d6SMark F. Adams    Options Database Key:
1117cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
1118ffc955d6SMark F. Adams 
1119ffc955d6SMark F. Adams    Level: intermediate
1120ffc955d6SMark F. Adams 
1121ffc955d6SMark F. Adams .seealso: ()
1122ffc955d6SMark F. Adams @*/
1123cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1124ffc955d6SMark F. Adams {
1125ffc955d6SMark F. Adams   PetscErrorCode ierr;
1126ffc955d6SMark F. Adams 
1127ffc955d6SMark F. Adams   PetscFunctionBegin;
1128ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1129cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1130ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1131ffc955d6SMark F. Adams }
1132ffc955d6SMark F. Adams 
1133cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1134ffc955d6SMark F. Adams {
1135ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1136ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1137ffc955d6SMark F. Adams 
1138ffc955d6SMark F. Adams   PetscFunctionBegin;
1139cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1140ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1141ffc955d6SMark F. Adams }
1142ffc955d6SMark F. Adams 
1143171cca9aSMark Adams /*@
1144cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1145171cca9aSMark Adams 
1146171cca9aSMark Adams    Collective on PC
1147171cca9aSMark Adams 
1148171cca9aSMark Adams    Input Parameters:
1149171cca9aSMark Adams +  pc - the preconditioner context
1150cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1151171cca9aSMark Adams 
1152171cca9aSMark Adams    Options Database Key:
1153cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
1154171cca9aSMark Adams 
1155171cca9aSMark Adams    Level: intermediate
1156171cca9aSMark Adams 
115739d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1158171cca9aSMark Adams @*/
1159171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1160171cca9aSMark Adams {
1161171cca9aSMark Adams   PetscErrorCode ierr;
1162171cca9aSMark Adams 
1163171cca9aSMark Adams   PetscFunctionBegin;
1164171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1165171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1166171cca9aSMark Adams   PetscFunctionReturn(0);
1167171cca9aSMark Adams }
1168171cca9aSMark Adams 
1169171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1170171cca9aSMark Adams {
1171171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1172171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1173171cca9aSMark Adams 
1174171cca9aSMark Adams   PetscFunctionBegin;
1175171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1176ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1177ffc955d6SMark F. Adams }
1178ffc955d6SMark F. Adams 
11794ef23d27SMark F. Adams /*@
1180ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1181ce7c7f2fSMark Adams 
1182ce7c7f2fSMark Adams    Collective on PC
1183ce7c7f2fSMark Adams 
1184ce7c7f2fSMark Adams    Input Parameters:
1185ce7c7f2fSMark Adams +  pc - the preconditioner context
1186ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1187ce7c7f2fSMark Adams 
1188ce7c7f2fSMark Adams    Options Database Key:
1189ce7c7f2fSMark Adams .  -pc_gamg_cpu_pin_coarse_grids
1190ce7c7f2fSMark Adams 
1191ce7c7f2fSMark Adams    Level: intermediate
1192ce7c7f2fSMark Adams 
119339d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1194ce7c7f2fSMark Adams @*/
1195ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1196ce7c7f2fSMark Adams {
1197ce7c7f2fSMark Adams   PetscErrorCode ierr;
1198ce7c7f2fSMark Adams 
1199ce7c7f2fSMark Adams   PetscFunctionBegin;
1200ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1201ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1202ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1203ce7c7f2fSMark Adams }
1204ce7c7f2fSMark Adams 
1205ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1206ce7c7f2fSMark Adams {
1207ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1208ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1209ce7c7f2fSMark Adams 
1210ce7c7f2fSMark Adams   PetscFunctionBegin;
1211ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1212ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1213ce7c7f2fSMark Adams }
1214ce7c7f2fSMark Adams 
1215ce7c7f2fSMark Adams /*@
1216ce7c7f2fSMark Adams    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1217ce7c7f2fSMark Adams 
1218ce7c7f2fSMark Adams    Collective on PC
1219ce7c7f2fSMark Adams 
1220ce7c7f2fSMark Adams    Input Parameters:
1221ce7c7f2fSMark Adams +  pc - the preconditioner context
1222ce7c7f2fSMark Adams -  flg - Layout type
1223ce7c7f2fSMark Adams 
1224ce7c7f2fSMark Adams    Options Database Key:
1225ce7c7f2fSMark Adams .  -pc_gamg_coarse_grid_layout_type
1226ce7c7f2fSMark Adams 
1227ce7c7f2fSMark Adams    Level: intermediate
1228ce7c7f2fSMark Adams 
122939d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1230ce7c7f2fSMark Adams @*/
1231ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1232ce7c7f2fSMark Adams {
1233ce7c7f2fSMark Adams   PetscErrorCode ierr;
1234ce7c7f2fSMark Adams 
1235ce7c7f2fSMark Adams   PetscFunctionBegin;
1236ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1237ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1238ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1239ce7c7f2fSMark Adams }
1240ce7c7f2fSMark Adams 
1241ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1242ce7c7f2fSMark Adams {
1243ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1244ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1245ce7c7f2fSMark Adams 
1246ce7c7f2fSMark Adams   PetscFunctionBegin;
1247ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1248ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1249ce7c7f2fSMark Adams }
1250ce7c7f2fSMark Adams 
1251ce7c7f2fSMark Adams /*@
12521cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
12534ef23d27SMark F. Adams 
12544ef23d27SMark F. Adams    Not collective on PC
12554ef23d27SMark F. Adams 
12564ef23d27SMark F. Adams    Input Parameters:
12571cc46a46SBarry Smith +  pc - the preconditioner
12581cc46a46SBarry Smith -  n - the maximum number of levels to use
12594ef23d27SMark F. Adams 
12604ef23d27SMark F. Adams    Options Database Key:
12614ef23d27SMark F. Adams .  -pc_mg_levels
12624ef23d27SMark F. Adams 
12634ef23d27SMark F. Adams    Level: intermediate
12644ef23d27SMark F. Adams 
12654ef23d27SMark F. Adams .seealso: ()
12664ef23d27SMark F. Adams @*/
12674ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12684ef23d27SMark F. Adams {
12694ef23d27SMark F. Adams   PetscErrorCode ierr;
12704ef23d27SMark F. Adams 
12714ef23d27SMark F. Adams   PetscFunctionBegin;
12724ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12734ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12744ef23d27SMark F. Adams   PetscFunctionReturn(0);
12754ef23d27SMark F. Adams }
12764ef23d27SMark F. Adams 
12771e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12784ef23d27SMark F. Adams {
12794ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12804ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12814ef23d27SMark F. Adams 
12824ef23d27SMark F. Adams   PetscFunctionBegin;
12839d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12844ef23d27SMark F. Adams   PetscFunctionReturn(0);
12854ef23d27SMark F. Adams }
12864ef23d27SMark F. Adams 
12873542efc5SMark F. Adams /*@
12883542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12893542efc5SMark F. Adams 
12903542efc5SMark F. Adams    Not collective on PC
12913542efc5SMark F. Adams 
12923542efc5SMark F. Adams    Input Parameters:
12931cc46a46SBarry Smith +  pc - the preconditioner context
1294055c8bd0SJed Brown .  threshold - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
1295055c8bd0SJed Brown -  n - number of threshold values provided in array
12963542efc5SMark F. Adams 
12973542efc5SMark F. Adams    Options Database Key:
12981cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
12993542efc5SMark F. Adams 
130095452b02SPatrick Sanan    Notes:
1301af3c827dSMark 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.
1302af3c827dSMark 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.
1303cab9ed1eSBarry Smith 
1304055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1305055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1306055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1307055c8bd0SJed Brown 
13083542efc5SMark F. Adams    Level: intermediate
13093542efc5SMark F. Adams 
1310af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
13113542efc5SMark F. Adams @*/
1312c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
13133542efc5SMark F. Adams {
13143542efc5SMark F. Adams   PetscErrorCode ierr;
13153542efc5SMark F. Adams 
13163542efc5SMark F. Adams   PetscFunctionBegin;
13173542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1318055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
1319c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
13203542efc5SMark F. Adams   PetscFunctionReturn(0);
13213542efc5SMark F. Adams }
13223542efc5SMark F. Adams 
1323c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
13243542efc5SMark F. Adams {
1325c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1326c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1327c1eae691SMark Adams   PetscInt i;
1328c1eae691SMark Adams   PetscFunctionBegin;
1329055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1330055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1331c1eae691SMark Adams   PetscFunctionReturn(0);
1332c1eae691SMark Adams }
1333c1eae691SMark Adams 
1334c1eae691SMark Adams /*@
1335c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1336c1eae691SMark Adams 
1337c1eae691SMark Adams    Not collective on PC
1338c1eae691SMark Adams 
1339c1eae691SMark Adams    Input Parameters:
1340c1eae691SMark Adams +  pc - the preconditioner context
1341c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1342c1eae691SMark Adams 
1343c1eae691SMark Adams    Options Database Key:
1344c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1345c1eae691SMark Adams 
1346055c8bd0SJed Brown    Notes:
1347055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1348055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1349055c8bd0SJed Brown 
1350c1eae691SMark Adams    Level: advanced
1351c1eae691SMark Adams 
1352055c8bd0SJed Brown .seealso: PCGAMGSetThreshold()
1353c1eae691SMark Adams @*/
1354c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1355c1eae691SMark Adams {
1356c1eae691SMark Adams   PetscErrorCode ierr;
13573542efc5SMark F. Adams 
13583542efc5SMark F. Adams   PetscFunctionBegin;
1359c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1360c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1361c1eae691SMark Adams   PetscFunctionReturn(0);
1362c1eae691SMark Adams }
1363c1eae691SMark Adams 
1364c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1365c1eae691SMark Adams {
1366c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1367c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1368c1eae691SMark Adams   PetscFunctionBegin;
1369c1eae691SMark Adams   pc_gamg->threshold_scale = v;
13703542efc5SMark F. Adams   PetscFunctionReturn(0);
13713542efc5SMark F. Adams }
13723542efc5SMark F. Adams 
1373e20c40e8SBarry Smith /*@C
1374c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1375676e1743SMark F. Adams 
1376676e1743SMark F. Adams    Collective on PC
1377676e1743SMark F. Adams 
1378676e1743SMark F. Adams    Input Parameters:
1379c60c7ad4SBarry Smith +  pc - the preconditioner context
1380c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1381676e1743SMark F. Adams 
1382676e1743SMark F. Adams    Options Database Key:
1383cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1384676e1743SMark F. Adams 
1385676e1743SMark F. Adams    Level: intermediate
1386676e1743SMark F. Adams 
1387cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1388676e1743SMark F. Adams @*/
138919fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1390676e1743SMark F. Adams {
1391676e1743SMark F. Adams   PetscErrorCode ierr;
1392676e1743SMark F. Adams 
1393676e1743SMark F. Adams   PetscFunctionBegin;
1394676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1395806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1396676e1743SMark F. Adams   PetscFunctionReturn(0);
1397676e1743SMark F. Adams }
1398676e1743SMark F. Adams 
1399e20c40e8SBarry Smith /*@C
1400c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1401c60c7ad4SBarry Smith 
1402c60c7ad4SBarry Smith    Collective on PC
1403c60c7ad4SBarry Smith 
1404c60c7ad4SBarry Smith    Input Parameter:
1405c60c7ad4SBarry Smith .  pc - the preconditioner context
1406c60c7ad4SBarry Smith 
1407c60c7ad4SBarry Smith    Output Parameter:
1408c60c7ad4SBarry Smith .  type - the type of algorithm used
1409c60c7ad4SBarry Smith 
1410c60c7ad4SBarry Smith    Level: intermediate
1411c60c7ad4SBarry Smith 
14121c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1413c60c7ad4SBarry Smith @*/
1414c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1415c60c7ad4SBarry Smith {
1416c60c7ad4SBarry Smith   PetscErrorCode ierr;
1417c60c7ad4SBarry Smith 
1418c60c7ad4SBarry Smith   PetscFunctionBegin;
1419c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1420c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1421c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1422c60c7ad4SBarry Smith }
1423c60c7ad4SBarry Smith 
1424c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1425c60c7ad4SBarry Smith {
1426c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1427c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1428c60c7ad4SBarry Smith 
1429c60c7ad4SBarry Smith   PetscFunctionBegin;
1430c60c7ad4SBarry Smith   *type = pc_gamg->type;
1431c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1432c60c7ad4SBarry Smith }
1433c60c7ad4SBarry Smith 
14341e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1435676e1743SMark F. Adams {
14369d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
14371ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
14381ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1439676e1743SMark F. Adams 
1440676e1743SMark F. Adams   PetscFunctionBegin;
1441c60c7ad4SBarry Smith   pc_gamg->type = type;
14421c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
14439d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
14441ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
14451ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
14461ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1447e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14483ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
14493ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
14503ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
14513ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
14523ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
14533ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
14543ae0bb68SMark Adams     pc_gamg->data_sz = 0;
14551ab5ffc9SJed Brown   }
14561ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
14571ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
14589d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1459676e1743SMark F. Adams   PetscFunctionReturn(0);
1460676e1743SMark F. Adams }
1461676e1743SMark F. Adams 
14629d766c59SMark Adams /* -------------------------------------------------------------------------- */
14639d766c59SMark Adams /*
14649d766c59SMark Adams    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
14659d766c59SMark Adams 
14669d766c59SMark Adams    Input Parameter:
14679d766c59SMark Adams .  pc - the preconditioner context
14689d766c59SMark Adams 
14699d766c59SMark Adams    Output Parameter:
14709d766c59SMark Adams .  gc - grid complexity = sum_i(nnz_i) / nnz_0
14719d766c59SMark Adams 
14729d766c59SMark Adams    Level: advanced
14739d766c59SMark Adams */
14749d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
14759d766c59SMark Adams {
14769d766c59SMark Adams   PetscErrorCode ierr;
14779d766c59SMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
14789d766c59SMark Adams   PC_MG_Levels   **mglevels = mg->levels;
14793966268fSBarry Smith   PetscInt       lev;
14803966268fSBarry Smith   PetscLogDouble nnz0 = 0, sgc = 0;
14819d766c59SMark Adams   MatInfo        info;
14823966268fSBarry Smith 
14839d766c59SMark Adams   PetscFunctionBegin;
1484dbf6bb8dSprj-   if (!pc->setupcalled) {
1485dbf6bb8dSprj-     *gc = 0;
1486dbf6bb8dSprj-     PetscFunctionReturn(0);
1487dbf6bb8dSprj-   }
14889d766c59SMark Adams   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
14893966268fSBarry Smith   for (lev=0; lev<mg->nlevels; lev++) {
149062a6e064SMark Adams     Mat dB;
149162a6e064SMark Adams     ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr);
149262a6e064SMark Adams     ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
14933966268fSBarry Smith     sgc += info.nz_used;
14949d766c59SMark Adams     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
14959d766c59SMark Adams   }
14963966268fSBarry Smith   if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0);
14973966268fSBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
14989d766c59SMark Adams   PetscFunctionReturn(0);
14999d766c59SMark Adams }
15009d766c59SMark Adams 
15015adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
15025adeb434SBarry Smith {
1503c1eae691SMark Adams   PetscErrorCode ierr,i;
15045adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
15055adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
150623b2d91dSMark Adams   PetscReal       gc=0;
15075adeb434SBarry Smith   PetscFunctionBegin;
15085adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1509459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1510c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1511c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1512c1eae691SMark Adams   }
1513459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1514459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1515cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1516cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1517cab9ed1eSBarry Smith   }
1518171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1519171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1520171cca9aSMark Adams   }
1521ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1522ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1523ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1524ce7c7f2fSMark Adams   }
1525ce7c7f2fSMark Adams #endif
1526ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1527ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1528ce7c7f2fSMark Adams   /* } else { */
1529ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1530ce7c7f2fSMark Adams   /* } */
15315adeb434SBarry Smith   if (pc_gamg->ops->view) {
15325adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
15335adeb434SBarry Smith   }
15349d766c59SMark Adams   ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr);
15359d766c59SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);CHKERRQ(ierr);
15365adeb434SBarry Smith   PetscFunctionReturn(0);
15375adeb434SBarry Smith }
15385adeb434SBarry Smith 
15394416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
15405b89ad90SMark F. Adams {
1541676e1743SMark F. Adams   PetscErrorCode ierr;
1542676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1543676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
154418c3aa7eSMark   PetscBool      flag,f2;
15453b4367a7SBarry Smith   MPI_Comm       comm;
154618c3aa7eSMark   char           prefix[256],tname[32];
1547c1eae691SMark Adams   PetscInt       i,n;
154814a9496bSBarry Smith   const char     *pcpre;
15490a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
15505b89ad90SMark F. Adams   PetscFunctionBegin;
15513b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1552e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
15531a1c1e04SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1554bd94a7aaSJed Brown   if (flag) {
1555bd94a7aaSJed Brown     ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
15561ab5ffc9SJed Brown   }
155718c3aa7eSMark   ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr);
155818c3aa7eSMark   if (flag) {
155918c3aa7eSMark     ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr);
156018c3aa7eSMark   }
1561cab9ed1eSBarry Smith   ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
156218c3aa7eSMark   f2 = PETSC_TRUE;
156318c3aa7eSMark   ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr);
156418c3aa7eSMark   if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0;
15651cc46a46SBarry Smith   ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1566a303c832SJed 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);
1567cf8ae1d3SMark 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);
1568ce7c7f2fSMark 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);
1569a0095786SMark   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);
157094ae4db5SBarry 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);
157118c3aa7eSMark   ierr = PetscOptionsInt("-pc_gamg_esteig_ksp_max_it","Number of iterations of eigen estimator","PCGAMGSetEstEigKSPMaxIt",pc_gamg->esteig_max_it,&pc_gamg->esteig_max_it,NULL);CHKERRQ(ierr);
157294ae4db5SBarry 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);
1573a303c832SJed 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);
157418c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
1575c1eae691SMark Adams   ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
157618c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1577efd3c5ceSMark Adams     if (!flag) n = 1;
1578c1eae691SMark Adams     i = n;
157918c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1580c1eae691SMark Adams   }
158194ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
158218c3aa7eSMark   {
158318c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
158418c3aa7eSMark     n = 2;
158518c3aa7eSMark     ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr);
158618c3aa7eSMark     if (flag) {
158718c3aa7eSMark       if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
158818c3aa7eSMark       ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr);
158918c3aa7eSMark     }
159018c3aa7eSMark   }
1591b7cbab4eSMark Adams   /* set options for subtype */
1592e55864a3SBarry Smith   if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
159318c3aa7eSMark 
159414a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
159514a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1596676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
15975b89ad90SMark F. Adams   PetscFunctionReturn(0);
15985b89ad90SMark F. Adams }
15995b89ad90SMark F. Adams 
16005b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
16015b89ad90SMark F. Adams /*MC
16021cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
16035b89ad90SMark F. Adams 
1604280d9858SJed Brown    Options Database Keys:
1605cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1606cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1607cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1608cab9ed1eSBarry 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
1609cab9ed1eSBarry 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>
1610cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1611cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
16126008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1613c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1614cab9ed1eSBarry Smith 
1615cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1616cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1617cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1618cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1619cab9ed1eSBarry Smith 
1620db9745e2SBarry Smith    Multigrid options:
1621db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1622db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1623db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1624cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
16255b89ad90SMark F. Adams 
16261cc46a46SBarry Smith 
162795452b02SPatrick Sanan   Notes:
162895452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1629db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1630db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1631db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
16321cc46a46SBarry Smith 
16335b89ad90SMark F. Adams   Level: intermediate
1634280d9858SJed Brown 
16351cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
163618c3aa7eSMark            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType()
16375b89ad90SMark F. Adams M*/
1638b2573a8aSBarry Smith 
16398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
16405b89ad90SMark F. Adams {
1641c1eae691SMark Adams   PetscErrorCode ierr,i;
16425b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
16435b89ad90SMark F. Adams   PC_MG          *mg;
16445b89ad90SMark F. Adams 
16455b89ad90SMark F. Adams   PetscFunctionBegin;
16461c1aac46SBarry Smith    /* register AMG type */
16471c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
16481c1aac46SBarry Smith 
16495b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
16501c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
16515b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
16525b89ad90SMark F. Adams 
16535b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1654b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
165569aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
16565b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
16575b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
16585b89ad90SMark F. Adams 
1659b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
16601ab5ffc9SJed Brown 
16619d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
16629d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
16639d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
16640a545947SLisandro Dalcin   pc_gamg->data    = NULL;
16655b89ad90SMark F. Adams 
16669d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
16675b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
16685b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
16695b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
16705b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
16715adeb434SBarry Smith   mg->view                = PCView_GAMG;
16725b89ad90SMark F. Adams 
167397d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
167497d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1675bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1676bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1677cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
167818c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr);
167918c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr);
168018c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr);
168118c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr);
16821cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1683cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1684171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1685ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1686ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1687bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1688c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1689bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1690c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1691bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
16929d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1693d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
16940c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1695171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1696a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1697a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1698038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
169925a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
170018c3aa7eSMark   for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1701c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
170218c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
17039ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1704d24ecf33SMark   ierr = PetscStrcpy(pc_gamg->esteig_type,NULL);CHKERRQ(ierr);
170518c3aa7eSMark   pc_gamg->esteig_max_it    = 10;
170618c3aa7eSMark   pc_gamg->use_sa_esteig    = -1;
170718c3aa7eSMark   pc_gamg->emin             = 0;
170818c3aa7eSMark   pc_gamg->emax             = 0;
170918c3aa7eSMark 
1710c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
17119d5b6da9SMark F. Adams 
1712bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1713bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
17145b89ad90SMark F. Adams   PetscFunctionReturn(0);
17155b89ad90SMark F. Adams }
17163e3471ccSMark Adams 
17173e3471ccSMark Adams /*@C
17183e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
17198a690491SBarry Smith     from PCInitializePackage().
17203e3471ccSMark Adams 
17213e3471ccSMark Adams  Level: developer
17223e3471ccSMark Adams 
17233e3471ccSMark Adams  .seealso: PetscInitialize()
17243e3471ccSMark Adams @*/
17253e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
17263e3471ccSMark Adams {
17273e3471ccSMark Adams   PetscErrorCode ierr;
17283e3471ccSMark Adams 
17293e3471ccSMark Adams   PetscFunctionBegin;
17303e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
17313e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
17323e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
17333e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
17348e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
17353e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1736c1c463dbSMark Adams 
1737c1c463dbSMark Adams   /* general events */
1738fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1739fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1740fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1741fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1742c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1743c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1744fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1745c1c463dbSMark Adams 
17465b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
17475b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
17485b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
17495b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
17505b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
17515b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
17525b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
17535b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
17545b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1755bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
17565b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
17575b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
17585b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
17595b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
17605b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
17615b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
17625b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
17635b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
17645b89ad90SMark F. Adams 
17655b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
17665b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
17675b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
17685b89ad90SMark F. Adams   /* create timer stages */
17695b89ad90SMark F. Adams #if defined GAMG_STAGES
17705b89ad90SMark F. Adams   {
17715b89ad90SMark F. Adams     char     str[32];
17725b89ad90SMark F. Adams     PetscInt lidx;
17735b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
17745b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
17755b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
17765b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
17775b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
17785b89ad90SMark F. Adams     }
17795b89ad90SMark F. Adams   }
17805b89ad90SMark F. Adams #endif
17815b89ad90SMark F. Adams #endif
17823e3471ccSMark Adams   PetscFunctionReturn(0);
17833e3471ccSMark Adams }
17843e3471ccSMark Adams 
17853e3471ccSMark Adams /*@C
17861c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
17871c1aac46SBarry Smith     called from PetscFinalize() automatically.
17883e3471ccSMark Adams 
17893e3471ccSMark Adams  Level: developer
17903e3471ccSMark Adams 
17913e3471ccSMark Adams  .seealso: PetscFinalize()
17923e3471ccSMark Adams @*/
17933e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
17943e3471ccSMark Adams {
17953e3471ccSMark Adams   PetscErrorCode ierr;
17963e3471ccSMark Adams 
17973e3471ccSMark Adams   PetscFunctionBegin;
17983e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
17993e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
18003e3471ccSMark Adams   PetscFunctionReturn(0);
18013e3471ccSMark Adams }
1802a36cf38bSToby Isaac 
1803a36cf38bSToby Isaac /*@C
1804a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1805a36cf38bSToby Isaac 
1806a36cf38bSToby Isaac  Input Parameters:
1807a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1808a36cf38bSToby Isaac  - create - function for creating the gamg context.
1809a36cf38bSToby Isaac 
1810a36cf38bSToby Isaac   Level: advanced
1811a36cf38bSToby Isaac 
18121c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1813a36cf38bSToby Isaac @*/
1814a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1815a36cf38bSToby Isaac {
1816a36cf38bSToby Isaac   PetscErrorCode ierr;
1817a36cf38bSToby Isaac 
1818a36cf38bSToby Isaac   PetscFunctionBegin;
1819a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1820a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1821a36cf38bSToby Isaac   PetscFunctionReturn(0);
1822a36cf38bSToby Isaac }
1823a36cf38bSToby Isaac 
1824