xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 0ed2132d46fdf59546ff8da87c5a9fab18cd3da0)
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 
4284b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
4294b1575e2SStefano Zampini {
4304b1575e2SStefano Zampini   PetscErrorCode ierr;
4314b1575e2SStefano Zampini   const char     *prefix;
4324b1575e2SStefano Zampini   char           addp[32];
4334b1575e2SStefano Zampini   PC_MG          *mg      = (PC_MG*)a_pc->data;
4344b1575e2SStefano Zampini   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4354b1575e2SStefano Zampini 
4364b1575e2SStefano Zampini   PetscFunctionBegin;
4374b1575e2SStefano Zampini   ierr = PCGetOptionsPrefix(a_pc,&prefix);CHKERRQ(ierr);
4384b1575e2SStefano Zampini   ierr = PetscInfo1(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);CHKERRQ(ierr);
4394b1575e2SStefano Zampini   ierr = MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);CHKERRQ(ierr);
4404b1575e2SStefano Zampini   ierr = MatSetOptionsPrefix(*Gmat2,prefix);CHKERRQ(ierr);
4414b1575e2SStefano Zampini   ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);CHKERRQ(ierr);
4424b1575e2SStefano Zampini   ierr = MatAppendOptionsPrefix(*Gmat2,addp);CHKERRQ(ierr);
4434b1575e2SStefano Zampini   ierr = MatProductSetType(*Gmat2,MATPRODUCT_AtB);CHKERRQ(ierr);
4444b1575e2SStefano Zampini   ierr = MatProductSetFromOptions(*Gmat2);CHKERRQ(ierr);
4454b1575e2SStefano Zampini   ierr = MatProductSymbolic(*Gmat2);CHKERRQ(ierr);
4464b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
4474b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
4484b1575e2SStefano Zampini   PetscFunctionReturn(0);
4494b1575e2SStefano Zampini }
4504b1575e2SStefano Zampini 
4515b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
4525b89ad90SMark F. Adams /*
4535b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4545b89ad90SMark F. Adams                     by setting data structures and options.
4555b89ad90SMark F. Adams 
4565b89ad90SMark F. Adams    Input Parameter:
4575b89ad90SMark F. Adams .  pc - the preconditioner context
4585b89ad90SMark F. Adams 
4595b89ad90SMark F. Adams */
4609d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
4615b89ad90SMark F. Adams {
4625b89ad90SMark F. Adams   PetscErrorCode ierr;
4639d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
4645b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4652adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
46618c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
4673b4367a7SBarry Smith   MPI_Comm       comm;
468c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
46918c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
47018c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
471a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
472569f4572SMark Adams   MatInfo        info;
473171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
4745ef31b24SMark F. Adams 
4755b89ad90SMark F. Adams   PetscFunctionBegin;
4763b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4773b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4783b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
479dfd5c07aSMark F. Adams 
48084d3f75bSMark F. Adams   if (pc_gamg->setup_count++ > 0) {
4811c1aac46SBarry Smith     if ((PetscBool)(!pc_gamg->reuse_prol)) {
482878e152fSMark F. Adams       /* reset everything */
483878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
484878e152fSMark F. Adams       pc->setupcalled = 0;
485806fa848SBarry Smith     } else {
48684d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
48703a628feSMark F. Adams       /* just do Galerkin grids */
48858471d46SMark F. Adams       Mat          B,dA,dB;
48958471d46SMark F. Adams 
49071959b99SBarry Smith       if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
4919d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
49258471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
49323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
49458471d46SMark F. Adams         /* (re)set to get dirty flag */
49523ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
49658471d46SMark F. Adams 
4972fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
498ef3f0257SMark Adams           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
499ef3f0257SMark Adams           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
500ef3f0257SMark Adams             ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
501ef3f0257SMark Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr);
502084a8fe3SJed Brown             ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);
50303a628feSMark F. Adams             mglevels[level]->A = B;
504806fa848SBarry Smith           } else {
505ef3f0257SMark Adams             ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr);
50623ee1639SBarry Smith             ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
50758471d46SMark F. Adams             ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
50803a628feSMark F. Adams           }
50923ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
51058471d46SMark F. Adams           dB   = B;
51158471d46SMark F. Adams         }
5125f8cf99dSMark F. Adams       }
513d5280255SMark F. Adams 
514d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
51558471d46SMark F. Adams       PetscFunctionReturn(0);
516eb07cef2SMark F. Adams     }
517878e152fSMark F. Adams   }
518f6536408SMark F. Adams 
519878e152fSMark F. Adams   if (!pc_gamg->data) {
520878e152fSMark F. Adams     if (pc_gamg->orig_data) {
521878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5220298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5232fa5cd67SKarl Rupp 
524878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
525878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
526878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5272fa5cd67SKarl Rupp 
528785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
529878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
530806fa848SBarry Smith     } else {
5311ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5327700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5339d5b6da9SMark F. Adams     }
534878e152fSMark F. Adams   }
535878e152fSMark F. Adams 
536878e152fSMark F. Adams   /* cache original data for reuse */
5371c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
538785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
539878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
540878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
541878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
542878e152fSMark F. Adams   }
543038e3b61SMark F. Adams 
544302f38e8SMark F. Adams   /* get basic dims */
545302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
546171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
54784d3f75bSMark F. Adams 
548569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
549569f4572SMark Adams   nnz0   = info.nz_used;
550569f4572SMark Adams   nnztot = info.nz_used;
55162294041SBarry 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);
552569f4572SMark Adams 
553a2f3521dSMark F. Adams   /* Get A_i and R_i */
55462294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
5559ab59c8bSMark Adams     pc_gamg->current_level = level;
55618c3aa7eSMark     if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
5575b89ad90SMark F. Adams     level1 = level + 1;
5580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
5590cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
560a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
561a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
562b4fbaa2aSMark F. Adams #endif
563a2f3521dSMark F. Adams #endif
564c8b0795cSMark F. Adams     { /* construct prolongator */
565725b86d8SJed Brown       Mat              Gmat;
5660cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
5677700e67bSMark Adams       Mat              Prol11;
568c8b0795cSMark F. Adams 
5697700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
5701ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
5717700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
572c8b0795cSMark F. Adams 
573a2f3521dSMark F. Adams       /* could have failed to create new level */
574a2f3521dSMark F. Adams       if (Prol11) {
575f7df55f0SStefano Zampini         const char *prefix;
576f7df55f0SStefano Zampini         char       addp[32];
577f7df55f0SStefano Zampini 
5789d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
5790298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
580a2f3521dSMark F. Adams 
581fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
582c8b0795cSMark F. Adams           /* smooth */
583fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
584c8b0795cSMark F. Adams         }
585c8b0795cSMark F. Adams 
5860c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
5871b18a24aSMark Adams           PetscInt bs;
5881b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
5890a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
590ffc955d6SMark F. Adams         }
591ffc955d6SMark F. Adams 
592f7df55f0SStefano Zampini         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
593f7df55f0SStefano Zampini         ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr);
594f7df55f0SStefano Zampini         ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",level);CHKERRQ(ierr);
595f7df55f0SStefano Zampini         ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr);
596f7df55f0SStefano Zampini         /* if we reuse the prolongator, then it is better to generate the transpose by default with CUDA
597f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
598f7df55f0SStefano Zampini #if defined(PETSC_HAVE_CUDA)
599f7df55f0SStefano Zampini         {
600f7df55f0SStefano Zampini           PetscBool ismpiaij;
601f7df55f0SStefano Zampini 
602f7df55f0SStefano Zampini           ierr = PetscObjectBaseTypeCompare((PetscObject)Prol11,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
603f7df55f0SStefano Zampini           if (ismpiaij) {
604f7df55f0SStefano Zampini             Mat Prol_d,Prol_o;
605f7df55f0SStefano Zampini 
606f7df55f0SStefano Zampini             ierr = MatMPIAIJGetSeqAIJ(Prol11,&Prol_d,&Prol_o,NULL);CHKERRQ(ierr);
607f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol_d,pc_gamg->reuse_prol);CHKERRQ(ierr);
608f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol_o,pc_gamg->reuse_prol);CHKERRQ(ierr);
609f7df55f0SStefano Zampini           } else {
610f7df55f0SStefano Zampini             ierr = MatSeqAIJCUSPARSESetGenerateTranspose(Prol11,pc_gamg->reuse_prol);CHKERRQ(ierr);
611f7df55f0SStefano Zampini           }
612f7df55f0SStefano Zampini         }
613f7df55f0SStefano Zampini #endif
614f7df55f0SStefano Zampini         ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr);
6154bde40a0SMark Adams         Parr[level1] = Prol11;
6164bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6174bde40a0SMark Adams 
618a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
61941b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
620a2f3521dSMark F. Adams     } /* construct prolongator scope */
6210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6220cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
623c8b0795cSMark F. Adams #endif
6247f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
625171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
626569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
62762294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
628a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
629a90e85d9SMark Adams #endif
630c8b0795cSMark F. Adams       break;
631c8b0795cSMark F. Adams     }
6320cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6330cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
634b4fbaa2aSMark F. Adams #endif
635171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
636171cca9aSMark Adams     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
637171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6380e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
639171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
640a2f3521dSMark F. Adams 
6410cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6420cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
643b4fbaa2aSMark F. Adams #endif
644171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
645569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
646569f4572SMark Adams     nnztot += info.nz_used;
6471d5b2942SMark 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);
648569f4572SMark Adams 
6490cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
650b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
651b4fbaa2aSMark F. Adams #endif
652a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6539ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
6549ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
655a90e85d9SMark Adams       level++;
656a90e85d9SMark Adams       break;
657a90e85d9SMark Adams     }
658c8b0795cSMark F. Adams   } /* levels */
659c8b0795cSMark F. Adams   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
660c8b0795cSMark F. Adams 
661569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
6629d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6635b89ad90SMark F. Adams   fine_level       = level;
6640298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
6655b89ad90SMark F. Adams 
66662294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
667*0ed2132dSStefano Zampini     PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP);
668*0ed2132dSStefano Zampini 
669d5280255SMark F. Adams     /* set default smoothers & set operators */
67062294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
671ffc955d6SMark F. Adams       KSP smoother;
672ffc955d6SMark F. Adams       PC  subpc;
673a2f3521dSMark F. Adams 
6749d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
675f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
676ffc955d6SMark F. Adams 
677a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
678a2f3521dSMark F. Adams       /* set ops */
67923ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
680a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
681a2f3521dSMark F. Adams 
682a2f3521dSMark F. Adams       /* set defaults */
6836c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
684a2f3521dSMark F. Adams 
6850c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
6860c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
6872d3561bbSSatish Balay         PetscInt sz;
6887a28f3e5SMark Adams         IS       *iss;
689a2f3521dSMark F. Adams 
6902d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
6917a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
6920c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
6930a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
6940c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
6957f66b68fSMark Adams         if (!sz) {
696ffc955d6SMark F. Adams           IS       is;
6970a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
6987a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
699a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
700806fa848SBarry Smith         } else {
701a94c3b12SMark F. Adams           PetscInt kk;
7027a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
703a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
7047a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
705a94c3b12SMark F. Adams           }
7067a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
707ffc955d6SMark F. Adams         }
7080298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
709ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
710806fa848SBarry Smith       } else {
711890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
712ffc955d6SMark F. Adams       }
713d5280255SMark F. Adams     }
714d5280255SMark F. Adams     {
715d5280255SMark F. Adams       /* coarse grid */
716d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
717d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
718*0ed2132dSStefano Zampini 
719d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
72023ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
721cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
722d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
723d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
724d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
725d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
72671959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
72771959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
728d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
729d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7309dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7312fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
73208e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
7335b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7345b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7355b42dca8SJed Brown          * view every subdomain as though they were different. */
7365b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
737d5280255SMark F. Adams       }
738cf8ae1d3SMark Adams     }
739d5280255SMark F. Adams 
740d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
741d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
742e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
743d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
74469aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
745d5280255SMark F. Adams 
74618c3aa7eSMark     /* setup cheby eigen estimates from SA */
747*0ed2132dSStefano Zampini     if (pc_gamg->use_sa_esteig==1) {
74818c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
74918c3aa7eSMark         KSP       smoother;
75018c3aa7eSMark         PetscBool ischeb;
751*0ed2132dSStefano Zampini 
752*0ed2132dSStefano Zampini         savesetfromoptions[level] = NULL;
75318c3aa7eSMark         ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
75418c3aa7eSMark         ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr);
75518c3aa7eSMark         if (ischeb) {
75618c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
757*0ed2132dSStefano Zampini 
758*0ed2132dSStefano Zampini           ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr); /* let command line emax override using SA's eigenvalues */
759*0ed2132dSStefano Zampini           if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) {
76018c3aa7eSMark             PC        subpc;
76118c3aa7eSMark             PetscBool isjac;
76218c3aa7eSMark             ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
76318c3aa7eSMark             ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr);
764*0ed2132dSStefano Zampini             if (isjac && pc_gamg->use_sa_esteig==1) {
76518c3aa7eSMark               PetscReal emax,emin;
766*0ed2132dSStefano Zampini 
76718c3aa7eSMark               emin = mg->min_eigen_DinvA[level];
76818c3aa7eSMark               emax = mg->max_eigen_DinvA[level];
76918c3aa7eSMark               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);
77018c3aa7eSMark               cheb->emin_computed = emin;
77118c3aa7eSMark               cheb->emax_computed = emax;
77218c3aa7eSMark               ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr);
773*0ed2132dSStefano Zampini 
774*0ed2132dSStefano Zampini               /* We have set the eigenvalues and consumed the transformation values
775*0ed2132dSStefano Zampini                  prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG
776*0ed2132dSStefano Zampini                  below when setfromoptions will be called again */
777*0ed2132dSStefano Zampini               savesetfromoptions[level] = smoother->ops->setfromoptions;
778*0ed2132dSStefano Zampini               smoother->ops->setfromoptions = NULL;
77918c3aa7eSMark             }
78018c3aa7eSMark           }
78118c3aa7eSMark         }
78218c3aa7eSMark       }
783*0ed2132dSStefano Zampini     }
784*0ed2132dSStefano Zampini 
785*0ed2132dSStefano Zampini     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
786*0ed2132dSStefano Zampini 
787*0ed2132dSStefano Zampini     /* restore Chebyshev smoother for next calls */
788*0ed2132dSStefano Zampini     if (pc_gamg->use_sa_esteig==1) {
789*0ed2132dSStefano Zampini       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
790*0ed2132dSStefano Zampini         if (savesetfromoptions[level]) {
791*0ed2132dSStefano Zampini           KSP smoother;
792*0ed2132dSStefano Zampini           ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
793*0ed2132dSStefano Zampini           smoother->ops->setfromoptions = savesetfromoptions[level];
794*0ed2132dSStefano Zampini         }
795*0ed2132dSStefano Zampini       }
796*0ed2132dSStefano Zampini     }
79718c3aa7eSMark 
798d5280255SMark F. Adams     /* clean up */
799d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
800587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
801587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8025b89ad90SMark F. Adams     }
803806fa848SBarry Smith   } else {
8045f8cf99dSMark F. Adams     KSP smoother;
805*0ed2132dSStefano Zampini 
806302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
8079d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
80823ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8095f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8109d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8115f8cf99dSMark F. Adams   }
8125b89ad90SMark F. Adams   PetscFunctionReturn(0);
8135b89ad90SMark F. Adams }
8145b89ad90SMark F. Adams 
815eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8165b89ad90SMark F. Adams /*
8175b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8185b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8195b89ad90SMark F. Adams 
8205b89ad90SMark F. Adams    Input Parameter:
8215b89ad90SMark F. Adams .  pc - the preconditioner context
8225b89ad90SMark F. Adams 
8235b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8245b89ad90SMark F. Adams */
8255b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8265b89ad90SMark F. Adams {
8275b89ad90SMark F. Adams   PetscErrorCode ierr;
8285b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8295b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8305b89ad90SMark F. Adams 
8315b89ad90SMark F. Adams   PetscFunctionBegin;
8325b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8339b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8349b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
8359b8ffb57SJed Brown   }
8361ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
8371ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
8385b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8395b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8405b89ad90SMark F. Adams   PetscFunctionReturn(0);
8415b89ad90SMark F. Adams }
8425b89ad90SMark F. Adams 
843676e1743SMark F. Adams /*@
844cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
845676e1743SMark F. Adams 
8461cc46a46SBarry Smith    Logically Collective on PC
847676e1743SMark F. Adams 
848676e1743SMark F. Adams    Input Parameters:
8491cc46a46SBarry Smith +  pc - the preconditioner context
8501cc46a46SBarry Smith -  n - the number of equations
851676e1743SMark F. Adams 
852676e1743SMark F. Adams 
853676e1743SMark F. Adams    Options Database Key:
8541cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
855676e1743SMark F. Adams 
85695452b02SPatrick Sanan    Notes:
85795452b02SPatrick 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
858cab9ed1eSBarry Smith           that has degrees of freedom
859cab9ed1eSBarry Smith 
860676e1743SMark F. Adams    Level: intermediate
861676e1743SMark F. Adams 
8621c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim()
863676e1743SMark F. Adams @*/
864676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
865676e1743SMark F. Adams {
866676e1743SMark F. Adams   PetscErrorCode ierr;
867676e1743SMark F. Adams 
868676e1743SMark F. Adams   PetscFunctionBegin;
869676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
871676e1743SMark F. Adams   PetscFunctionReturn(0);
872676e1743SMark F. Adams }
873676e1743SMark F. Adams 
8741e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
875676e1743SMark F. Adams {
876c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
877c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
878676e1743SMark F. Adams 
879676e1743SMark F. Adams   PetscFunctionBegin;
8809d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
881676e1743SMark F. Adams   PetscFunctionReturn(0);
882676e1743SMark F. Adams }
883676e1743SMark F. Adams 
884389730f3SMark F. Adams /*@
885cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
886389730f3SMark F. Adams 
887389730f3SMark F. Adams  Collective on PC
888389730f3SMark F. Adams 
889389730f3SMark F. Adams    Input Parameters:
8901cc46a46SBarry Smith +  pc - the preconditioner context
8911cc46a46SBarry Smith -  n - maximum number of equations to aim for
892389730f3SMark F. Adams 
893389730f3SMark F. Adams    Options Database Key:
8941cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
895389730f3SMark F. Adams 
89674329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
89774329af1SBarry Smith      has less than 1000 unknowns.
89874329af1SBarry Smith 
899389730f3SMark F. Adams    Level: intermediate
900389730f3SMark F. Adams 
9011c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim()
902389730f3SMark F. Adams @*/
903389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
904389730f3SMark F. Adams {
905389730f3SMark F. Adams   PetscErrorCode ierr;
906389730f3SMark F. Adams 
907389730f3SMark F. Adams   PetscFunctionBegin;
908389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
909389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
910389730f3SMark F. Adams   PetscFunctionReturn(0);
911389730f3SMark F. Adams }
912389730f3SMark F. Adams 
9131e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
914389730f3SMark F. Adams {
915389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
916389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
917389730f3SMark F. Adams 
918389730f3SMark F. Adams   PetscFunctionBegin;
9199d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
920389730f3SMark F. Adams   PetscFunctionReturn(0);
921389730f3SMark F. Adams }
922389730f3SMark F. Adams 
923676e1743SMark F. Adams /*@
924cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
925676e1743SMark F. Adams 
926676e1743SMark F. Adams    Collective on PC
927676e1743SMark F. Adams 
928676e1743SMark F. Adams    Input Parameters:
9291cc46a46SBarry Smith +  pc - the preconditioner context
9301cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
931676e1743SMark F. Adams 
932676e1743SMark F. Adams    Options Database Key:
9331cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
934676e1743SMark F. Adams 
93595452b02SPatrick Sanan    Notes:
93695452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
937cab9ed1eSBarry Smith 
938676e1743SMark F. Adams    Level: intermediate
939676e1743SMark F. Adams 
940676e1743SMark F. Adams .seealso: ()
941676e1743SMark F. Adams @*/
942cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
943676e1743SMark F. Adams {
944676e1743SMark F. Adams   PetscErrorCode ierr;
945676e1743SMark F. Adams 
946676e1743SMark F. Adams   PetscFunctionBegin;
947676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
948cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
949676e1743SMark F. Adams   PetscFunctionReturn(0);
950676e1743SMark F. Adams }
951676e1743SMark F. Adams 
952cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
953676e1743SMark F. Adams {
954c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
955c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
956676e1743SMark F. Adams 
957676e1743SMark F. Adams   PetscFunctionBegin;
9589d5b6da9SMark F. Adams   pc_gamg->repart = n;
959676e1743SMark F. Adams   PetscFunctionReturn(0);
960676e1743SMark F. Adams }
961676e1743SMark F. Adams 
962dfd5c07aSMark F. Adams /*@
96318c3aa7eSMark    PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby)
96418c3aa7eSMark 
96518c3aa7eSMark    Collective on PC
96618c3aa7eSMark 
96718c3aa7eSMark    Input Parameters:
96818c3aa7eSMark +  pc - the preconditioner context
96918c3aa7eSMark -  n - number of its
97018c3aa7eSMark 
97118c3aa7eSMark    Options Database Key:
97218c3aa7eSMark .  -pc_gamg_esteig_ksp_max_it <its>
97318c3aa7eSMark 
97418c3aa7eSMark    Notes:
97518c3aa7eSMark 
97618c3aa7eSMark    Level: intermediate
97718c3aa7eSMark 
97818c3aa7eSMark .seealso: ()
97918c3aa7eSMark @*/
98018c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n)
98118c3aa7eSMark {
98218c3aa7eSMark   PetscErrorCode ierr;
98318c3aa7eSMark 
98418c3aa7eSMark   PetscFunctionBegin;
98518c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
98618c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
98718c3aa7eSMark   PetscFunctionReturn(0);
98818c3aa7eSMark }
98918c3aa7eSMark 
99018c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n)
99118c3aa7eSMark {
99218c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
99318c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
99418c3aa7eSMark 
99518c3aa7eSMark   PetscFunctionBegin;
99618c3aa7eSMark   pc_gamg->esteig_max_it = n;
99718c3aa7eSMark   PetscFunctionReturn(0);
99818c3aa7eSMark }
99918c3aa7eSMark 
100018c3aa7eSMark /*@
100118c3aa7eSMark    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother
100218c3aa7eSMark 
100318c3aa7eSMark    Collective on PC
100418c3aa7eSMark 
100518c3aa7eSMark    Input Parameters:
100618c3aa7eSMark +  pc - the preconditioner context
100718c3aa7eSMark -  n - number of its
100818c3aa7eSMark 
100918c3aa7eSMark    Options Database Key:
101018c3aa7eSMark .  -pc_gamg_use_sa_esteig <true,false>
101118c3aa7eSMark 
101218c3aa7eSMark    Notes:
101318c3aa7eSMark 
101418c3aa7eSMark    Level: intermediate
101518c3aa7eSMark 
101618c3aa7eSMark .seealso: ()
101718c3aa7eSMark @*/
101818c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
101918c3aa7eSMark {
102018c3aa7eSMark   PetscErrorCode ierr;
102118c3aa7eSMark 
102218c3aa7eSMark   PetscFunctionBegin;
102318c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
102418c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
102518c3aa7eSMark   PetscFunctionReturn(0);
102618c3aa7eSMark }
102718c3aa7eSMark 
1028*0ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
102918c3aa7eSMark {
103018c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
103118c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
103218c3aa7eSMark 
103318c3aa7eSMark   PetscFunctionBegin;
103418c3aa7eSMark   pc_gamg->use_sa_esteig = n ? 1 : 0;
103518c3aa7eSMark   PetscFunctionReturn(0);
103618c3aa7eSMark }
103718c3aa7eSMark 
103818c3aa7eSMark /*@C
103918c3aa7eSMark    PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby)
104018c3aa7eSMark 
104118c3aa7eSMark    Collective on PC
104218c3aa7eSMark 
104318c3aa7eSMark    Input Parameters:
104418c3aa7eSMark +  pc - the preconditioner context
104518c3aa7eSMark -  t - ksp type
104618c3aa7eSMark 
104718c3aa7eSMark    Options Database Key:
104818c3aa7eSMark .  -pc_gamg_esteig_ksp_type <type>
104918c3aa7eSMark 
105018c3aa7eSMark    Notes:
105118c3aa7eSMark 
105218c3aa7eSMark    Level: intermediate
105318c3aa7eSMark 
105418c3aa7eSMark .seealso: ()
105518c3aa7eSMark @*/
105618c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[])
105718c3aa7eSMark {
105818c3aa7eSMark   PetscErrorCode ierr;
105918c3aa7eSMark 
106018c3aa7eSMark   PetscFunctionBegin;
106118c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
106218c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr);
106318c3aa7eSMark   PetscFunctionReturn(0);
106418c3aa7eSMark }
106518c3aa7eSMark 
106618c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[])
106718c3aa7eSMark {
106818c3aa7eSMark   PetscErrorCode ierr;
106918c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
107018c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
107118c3aa7eSMark 
107218c3aa7eSMark   PetscFunctionBegin;
107318c3aa7eSMark   ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr);
107418c3aa7eSMark   PetscFunctionReturn(0);
107518c3aa7eSMark }
107618c3aa7eSMark 
107718c3aa7eSMark /*@
107818c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
107918c3aa7eSMark 
108018c3aa7eSMark    Collective on PC
108118c3aa7eSMark 
108218c3aa7eSMark    Input Parameters:
108318c3aa7eSMark +  pc - the preconditioner context
108418c3aa7eSMark -  emax - max eigenvalue
108518c3aa7eSMark .  emin - min eigenvalue
108618c3aa7eSMark 
108718c3aa7eSMark    Options Database Key:
108818c3aa7eSMark .  -pc_gamg_eigenvalues
108918c3aa7eSMark 
109018c3aa7eSMark    Level: intermediate
109118c3aa7eSMark 
109218c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType()
109318c3aa7eSMark @*/
109418c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
109518c3aa7eSMark {
109618c3aa7eSMark   PetscErrorCode ierr;
109718c3aa7eSMark 
109818c3aa7eSMark   PetscFunctionBegin;
109918c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
110018c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr);
110118c3aa7eSMark   PetscFunctionReturn(0);
110218c3aa7eSMark }
110318c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
110418c3aa7eSMark {
110518c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
110618c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
110718c3aa7eSMark 
110818c3aa7eSMark   PetscFunctionBegin;
110918c3aa7eSMark   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);
111018c3aa7eSMark   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);
111118c3aa7eSMark   pc_gamg->emax = emax;
111218c3aa7eSMark   pc_gamg->emin = emin;
111318c3aa7eSMark 
111418c3aa7eSMark   PetscFunctionReturn(0);
111518c3aa7eSMark }
111618c3aa7eSMark 
111718c3aa7eSMark /*@
1118cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1119dfd5c07aSMark F. Adams 
1120dfd5c07aSMark F. Adams    Collective on PC
1121dfd5c07aSMark F. Adams 
1122dfd5c07aSMark F. Adams    Input Parameters:
11231cc46a46SBarry Smith +  pc - the preconditioner context
11241cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1125dfd5c07aSMark F. Adams 
1126dfd5c07aSMark F. Adams    Options Database Key:
11271cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1128dfd5c07aSMark F. Adams 
1129dfd5c07aSMark F. Adams    Level: intermediate
1130dfd5c07aSMark F. Adams 
113195452b02SPatrick Sanan    Notes:
113295452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1133cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
1134cab9ed1eSBarry Smith 
1135dfd5c07aSMark F. Adams .seealso: ()
1136dfd5c07aSMark F. Adams @*/
11371cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1138dfd5c07aSMark F. Adams {
1139dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1140dfd5c07aSMark F. Adams 
1141dfd5c07aSMark F. Adams   PetscFunctionBegin;
1142dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11431cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1144dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1145dfd5c07aSMark F. Adams }
1146dfd5c07aSMark F. Adams 
11471cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1148dfd5c07aSMark F. Adams {
1149dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1150dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1151dfd5c07aSMark F. Adams 
1152dfd5c07aSMark F. Adams   PetscFunctionBegin;
1153dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1154dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1155dfd5c07aSMark F. Adams }
1156dfd5c07aSMark F. Adams 
1157ffc955d6SMark F. Adams /*@
1158cab9ed1eSBarry 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.
1159ffc955d6SMark F. Adams 
1160ffc955d6SMark F. Adams    Collective on PC
1161ffc955d6SMark F. Adams 
1162ffc955d6SMark F. Adams    Input Parameters:
1163cab9ed1eSBarry Smith +  pc - the preconditioner context
1164cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1165ffc955d6SMark F. Adams 
1166ffc955d6SMark F. Adams    Options Database Key:
1167cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
1168ffc955d6SMark F. Adams 
1169ffc955d6SMark F. Adams    Level: intermediate
1170ffc955d6SMark F. Adams 
1171ffc955d6SMark F. Adams .seealso: ()
1172ffc955d6SMark F. Adams @*/
1173cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1174ffc955d6SMark F. Adams {
1175ffc955d6SMark F. Adams   PetscErrorCode ierr;
1176ffc955d6SMark F. Adams 
1177ffc955d6SMark F. Adams   PetscFunctionBegin;
1178ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1179cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1180ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1181ffc955d6SMark F. Adams }
1182ffc955d6SMark F. Adams 
1183cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1184ffc955d6SMark F. Adams {
1185ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1186ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1187ffc955d6SMark F. Adams 
1188ffc955d6SMark F. Adams   PetscFunctionBegin;
1189cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1190ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1191ffc955d6SMark F. Adams }
1192ffc955d6SMark F. Adams 
1193171cca9aSMark Adams /*@
1194cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1195171cca9aSMark Adams 
1196171cca9aSMark Adams    Collective on PC
1197171cca9aSMark Adams 
1198171cca9aSMark Adams    Input Parameters:
1199171cca9aSMark Adams +  pc - the preconditioner context
1200cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1201171cca9aSMark Adams 
1202171cca9aSMark Adams    Options Database Key:
1203cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
1204171cca9aSMark Adams 
1205171cca9aSMark Adams    Level: intermediate
1206171cca9aSMark Adams 
120739d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1208171cca9aSMark Adams @*/
1209171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1210171cca9aSMark Adams {
1211171cca9aSMark Adams   PetscErrorCode ierr;
1212171cca9aSMark Adams 
1213171cca9aSMark Adams   PetscFunctionBegin;
1214171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1215171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1216171cca9aSMark Adams   PetscFunctionReturn(0);
1217171cca9aSMark Adams }
1218171cca9aSMark Adams 
1219171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1220171cca9aSMark Adams {
1221171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1222171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1223171cca9aSMark Adams 
1224171cca9aSMark Adams   PetscFunctionBegin;
1225171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1226ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1227ffc955d6SMark F. Adams }
1228ffc955d6SMark F. Adams 
12294ef23d27SMark F. Adams /*@
1230ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1231ce7c7f2fSMark Adams 
1232ce7c7f2fSMark Adams    Collective on PC
1233ce7c7f2fSMark Adams 
1234ce7c7f2fSMark Adams    Input Parameters:
1235ce7c7f2fSMark Adams +  pc - the preconditioner context
1236ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1237ce7c7f2fSMark Adams 
1238ce7c7f2fSMark Adams    Options Database Key:
1239ce7c7f2fSMark Adams .  -pc_gamg_cpu_pin_coarse_grids
1240ce7c7f2fSMark Adams 
1241ce7c7f2fSMark Adams    Level: intermediate
1242ce7c7f2fSMark Adams 
124339d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1244ce7c7f2fSMark Adams @*/
1245ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1246ce7c7f2fSMark Adams {
1247ce7c7f2fSMark Adams   PetscErrorCode ierr;
1248ce7c7f2fSMark Adams 
1249ce7c7f2fSMark Adams   PetscFunctionBegin;
1250ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1251ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1252ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1253ce7c7f2fSMark Adams }
1254ce7c7f2fSMark Adams 
1255ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1256ce7c7f2fSMark Adams {
1257ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1258ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1259ce7c7f2fSMark Adams 
1260ce7c7f2fSMark Adams   PetscFunctionBegin;
1261ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1262ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1263ce7c7f2fSMark Adams }
1264ce7c7f2fSMark Adams 
1265ce7c7f2fSMark Adams /*@
1266ce7c7f2fSMark Adams    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1267ce7c7f2fSMark Adams 
1268ce7c7f2fSMark Adams    Collective on PC
1269ce7c7f2fSMark Adams 
1270ce7c7f2fSMark Adams    Input Parameters:
1271ce7c7f2fSMark Adams +  pc - the preconditioner context
1272ce7c7f2fSMark Adams -  flg - Layout type
1273ce7c7f2fSMark Adams 
1274ce7c7f2fSMark Adams    Options Database Key:
1275ce7c7f2fSMark Adams .  -pc_gamg_coarse_grid_layout_type
1276ce7c7f2fSMark Adams 
1277ce7c7f2fSMark Adams    Level: intermediate
1278ce7c7f2fSMark Adams 
127939d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1280ce7c7f2fSMark Adams @*/
1281ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1282ce7c7f2fSMark Adams {
1283ce7c7f2fSMark Adams   PetscErrorCode ierr;
1284ce7c7f2fSMark Adams 
1285ce7c7f2fSMark Adams   PetscFunctionBegin;
1286ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1287ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1288ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1289ce7c7f2fSMark Adams }
1290ce7c7f2fSMark Adams 
1291ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1292ce7c7f2fSMark Adams {
1293ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1294ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1295ce7c7f2fSMark Adams 
1296ce7c7f2fSMark Adams   PetscFunctionBegin;
1297ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1298ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1299ce7c7f2fSMark Adams }
1300ce7c7f2fSMark Adams 
1301ce7c7f2fSMark Adams /*@
13021cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
13034ef23d27SMark F. Adams 
13044ef23d27SMark F. Adams    Not collective on PC
13054ef23d27SMark F. Adams 
13064ef23d27SMark F. Adams    Input Parameters:
13071cc46a46SBarry Smith +  pc - the preconditioner
13081cc46a46SBarry Smith -  n - the maximum number of levels to use
13094ef23d27SMark F. Adams 
13104ef23d27SMark F. Adams    Options Database Key:
13114ef23d27SMark F. Adams .  -pc_mg_levels
13124ef23d27SMark F. Adams 
13134ef23d27SMark F. Adams    Level: intermediate
13144ef23d27SMark F. Adams 
13154ef23d27SMark F. Adams .seealso: ()
13164ef23d27SMark F. Adams @*/
13174ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
13184ef23d27SMark F. Adams {
13194ef23d27SMark F. Adams   PetscErrorCode ierr;
13204ef23d27SMark F. Adams 
13214ef23d27SMark F. Adams   PetscFunctionBegin;
13224ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13234ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
13244ef23d27SMark F. Adams   PetscFunctionReturn(0);
13254ef23d27SMark F. Adams }
13264ef23d27SMark F. Adams 
13271e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
13284ef23d27SMark F. Adams {
13294ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
13304ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13314ef23d27SMark F. Adams 
13324ef23d27SMark F. Adams   PetscFunctionBegin;
13339d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
13344ef23d27SMark F. Adams   PetscFunctionReturn(0);
13354ef23d27SMark F. Adams }
13364ef23d27SMark F. Adams 
13373542efc5SMark F. Adams /*@
13383542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
13393542efc5SMark F. Adams 
13403542efc5SMark F. Adams    Not collective on PC
13413542efc5SMark F. Adams 
13423542efc5SMark F. Adams    Input Parameters:
13431cc46a46SBarry Smith +  pc - the preconditioner context
1344055c8bd0SJed 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
1345055c8bd0SJed Brown -  n - number of threshold values provided in array
13463542efc5SMark F. Adams 
13473542efc5SMark F. Adams    Options Database Key:
13481cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
13493542efc5SMark F. Adams 
135095452b02SPatrick Sanan    Notes:
1351af3c827dSMark 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.
1352af3c827dSMark 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.
1353cab9ed1eSBarry Smith 
1354055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1355055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1356055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1357055c8bd0SJed Brown 
13583542efc5SMark F. Adams    Level: intermediate
13593542efc5SMark F. Adams 
1360af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
13613542efc5SMark F. Adams @*/
1362c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
13633542efc5SMark F. Adams {
13643542efc5SMark F. Adams   PetscErrorCode ierr;
13653542efc5SMark F. Adams 
13663542efc5SMark F. Adams   PetscFunctionBegin;
13673542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1368055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
1369c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
13703542efc5SMark F. Adams   PetscFunctionReturn(0);
13713542efc5SMark F. Adams }
13723542efc5SMark F. Adams 
1373c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
13743542efc5SMark F. Adams {
1375c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1376c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1377c1eae691SMark Adams   PetscInt i;
1378c1eae691SMark Adams   PetscFunctionBegin;
1379055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1380055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1381c1eae691SMark Adams   PetscFunctionReturn(0);
1382c1eae691SMark Adams }
1383c1eae691SMark Adams 
1384c1eae691SMark Adams /*@
1385c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1386c1eae691SMark Adams 
1387c1eae691SMark Adams    Not collective on PC
1388c1eae691SMark Adams 
1389c1eae691SMark Adams    Input Parameters:
1390c1eae691SMark Adams +  pc - the preconditioner context
1391c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1392c1eae691SMark Adams 
1393c1eae691SMark Adams    Options Database Key:
1394c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1395c1eae691SMark Adams 
1396055c8bd0SJed Brown    Notes:
1397055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1398055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1399055c8bd0SJed Brown 
1400c1eae691SMark Adams    Level: advanced
1401c1eae691SMark Adams 
1402055c8bd0SJed Brown .seealso: PCGAMGSetThreshold()
1403c1eae691SMark Adams @*/
1404c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1405c1eae691SMark Adams {
1406c1eae691SMark Adams   PetscErrorCode ierr;
14073542efc5SMark F. Adams 
14083542efc5SMark F. Adams   PetscFunctionBegin;
1409c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1410c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1411c1eae691SMark Adams   PetscFunctionReturn(0);
1412c1eae691SMark Adams }
1413c1eae691SMark Adams 
1414c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1415c1eae691SMark Adams {
1416c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1417c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1418c1eae691SMark Adams   PetscFunctionBegin;
1419c1eae691SMark Adams   pc_gamg->threshold_scale = v;
14203542efc5SMark F. Adams   PetscFunctionReturn(0);
14213542efc5SMark F. Adams }
14223542efc5SMark F. Adams 
1423e20c40e8SBarry Smith /*@C
1424c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1425676e1743SMark F. Adams 
1426676e1743SMark F. Adams    Collective on PC
1427676e1743SMark F. Adams 
1428676e1743SMark F. Adams    Input Parameters:
1429c60c7ad4SBarry Smith +  pc - the preconditioner context
1430c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1431676e1743SMark F. Adams 
1432676e1743SMark F. Adams    Options Database Key:
1433cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1434676e1743SMark F. Adams 
1435676e1743SMark F. Adams    Level: intermediate
1436676e1743SMark F. Adams 
1437cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1438676e1743SMark F. Adams @*/
143919fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1440676e1743SMark F. Adams {
1441676e1743SMark F. Adams   PetscErrorCode ierr;
1442676e1743SMark F. Adams 
1443676e1743SMark F. Adams   PetscFunctionBegin;
1444676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1445806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1446676e1743SMark F. Adams   PetscFunctionReturn(0);
1447676e1743SMark F. Adams }
1448676e1743SMark F. Adams 
1449e20c40e8SBarry Smith /*@C
1450c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1451c60c7ad4SBarry Smith 
1452c60c7ad4SBarry Smith    Collective on PC
1453c60c7ad4SBarry Smith 
1454c60c7ad4SBarry Smith    Input Parameter:
1455c60c7ad4SBarry Smith .  pc - the preconditioner context
1456c60c7ad4SBarry Smith 
1457c60c7ad4SBarry Smith    Output Parameter:
1458c60c7ad4SBarry Smith .  type - the type of algorithm used
1459c60c7ad4SBarry Smith 
1460c60c7ad4SBarry Smith    Level: intermediate
1461c60c7ad4SBarry Smith 
14621c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1463c60c7ad4SBarry Smith @*/
1464c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1465c60c7ad4SBarry Smith {
1466c60c7ad4SBarry Smith   PetscErrorCode ierr;
1467c60c7ad4SBarry Smith 
1468c60c7ad4SBarry Smith   PetscFunctionBegin;
1469c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1470c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1471c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1472c60c7ad4SBarry Smith }
1473c60c7ad4SBarry Smith 
1474c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1475c60c7ad4SBarry Smith {
1476c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1477c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1478c60c7ad4SBarry Smith 
1479c60c7ad4SBarry Smith   PetscFunctionBegin;
1480c60c7ad4SBarry Smith   *type = pc_gamg->type;
1481c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1482c60c7ad4SBarry Smith }
1483c60c7ad4SBarry Smith 
14841e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1485676e1743SMark F. Adams {
14869d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
14871ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
14881ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1489676e1743SMark F. Adams 
1490676e1743SMark F. Adams   PetscFunctionBegin;
1491c60c7ad4SBarry Smith   pc_gamg->type = type;
14921c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
14939d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
14941ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
14951ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
14961ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1497e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14983ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
14993ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
15003ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
15013ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
15023ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
15033ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
15043ae0bb68SMark Adams     pc_gamg->data_sz = 0;
15051ab5ffc9SJed Brown   }
15061ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
15071ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
15089d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1509676e1743SMark F. Adams   PetscFunctionReturn(0);
1510676e1743SMark F. Adams }
1511676e1743SMark F. Adams 
15129d766c59SMark Adams /* -------------------------------------------------------------------------- */
15139d766c59SMark Adams /*
15149d766c59SMark Adams    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
15159d766c59SMark Adams 
15169d766c59SMark Adams    Input Parameter:
15179d766c59SMark Adams .  pc - the preconditioner context
15189d766c59SMark Adams 
15199d766c59SMark Adams    Output Parameter:
15209d766c59SMark Adams .  gc - grid complexity = sum_i(nnz_i) / nnz_0
15219d766c59SMark Adams 
15229d766c59SMark Adams    Level: advanced
15239d766c59SMark Adams */
15249d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
15259d766c59SMark Adams {
15269d766c59SMark Adams   PetscErrorCode ierr;
15279d766c59SMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
15289d766c59SMark Adams   PC_MG_Levels   **mglevels = mg->levels;
15293966268fSBarry Smith   PetscInt       lev;
15303966268fSBarry Smith   PetscLogDouble nnz0 = 0, sgc = 0;
15319d766c59SMark Adams   MatInfo        info;
15323966268fSBarry Smith 
15339d766c59SMark Adams   PetscFunctionBegin;
1534dbf6bb8dSprj-   if (!pc->setupcalled) {
1535dbf6bb8dSprj-     *gc = 0;
1536dbf6bb8dSprj-     PetscFunctionReturn(0);
1537dbf6bb8dSprj-   }
15389d766c59SMark Adams   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
15393966268fSBarry Smith   for (lev=0; lev<mg->nlevels; lev++) {
154062a6e064SMark Adams     Mat dB;
154162a6e064SMark Adams     ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr);
154262a6e064SMark Adams     ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
15433966268fSBarry Smith     sgc += info.nz_used;
15449d766c59SMark Adams     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
15459d766c59SMark Adams   }
15463966268fSBarry Smith   if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0);
15473966268fSBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
15489d766c59SMark Adams   PetscFunctionReturn(0);
15499d766c59SMark Adams }
15509d766c59SMark Adams 
15515adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
15525adeb434SBarry Smith {
1553c1eae691SMark Adams   PetscErrorCode ierr,i;
15545adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
15555adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
155623b2d91dSMark Adams   PetscReal       gc=0;
15575adeb434SBarry Smith   PetscFunctionBegin;
15585adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1559459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1560c1eae691SMark Adams   for (i=0;i<pc_gamg->current_level;i++) {
1561c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1562c1eae691SMark Adams   }
1563459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1564459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1565cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1566cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1567cab9ed1eSBarry Smith   }
1568171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1569171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1570171cca9aSMark Adams   }
1571ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1572ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1573ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1574ce7c7f2fSMark Adams   }
1575ce7c7f2fSMark Adams #endif
1576ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1577ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1578ce7c7f2fSMark Adams   /* } else { */
1579ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1580ce7c7f2fSMark Adams   /* } */
15815adeb434SBarry Smith   if (pc_gamg->ops->view) {
15825adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
15835adeb434SBarry Smith   }
15849d766c59SMark Adams   ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr);
15859d766c59SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);CHKERRQ(ierr);
15865adeb434SBarry Smith   PetscFunctionReturn(0);
15875adeb434SBarry Smith }
15885adeb434SBarry Smith 
15894416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
15905b89ad90SMark F. Adams {
1591676e1743SMark F. Adams   PetscErrorCode ierr;
1592676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1593676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
159418c3aa7eSMark   PetscBool      flag,f2;
15953b4367a7SBarry Smith   MPI_Comm       comm;
159618c3aa7eSMark   char           prefix[256],tname[32];
1597c1eae691SMark Adams   PetscInt       i,n;
159814a9496bSBarry Smith   const char     *pcpre;
15990a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
16005b89ad90SMark F. Adams   PetscFunctionBegin;
16013b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1602e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
16031a1c1e04SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1604bd94a7aaSJed Brown   if (flag) {
1605bd94a7aaSJed Brown     ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
16061ab5ffc9SJed Brown   }
160718c3aa7eSMark   ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr);
160818c3aa7eSMark   if (flag) {
160918c3aa7eSMark     ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr);
161018c3aa7eSMark   }
1611cab9ed1eSBarry Smith   ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
161218c3aa7eSMark   f2 = PETSC_TRUE;
161318c3aa7eSMark   ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr);
161418c3aa7eSMark   if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0;
16151cc46a46SBarry Smith   ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1616a303c832SJed 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);
1617cf8ae1d3SMark 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);
1618ce7c7f2fSMark 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);
1619a0095786SMark   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);
162094ae4db5SBarry 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);
162118c3aa7eSMark   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);
162294ae4db5SBarry 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);
1623a303c832SJed 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);
162418c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
1625c1eae691SMark Adams   ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
162618c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1627efd3c5ceSMark Adams     if (!flag) n = 1;
1628c1eae691SMark Adams     i = n;
162918c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1630c1eae691SMark Adams   }
163194ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
163218c3aa7eSMark   {
163318c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
163418c3aa7eSMark     n = 2;
163518c3aa7eSMark     ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr);
163618c3aa7eSMark     if (flag) {
163718c3aa7eSMark       if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
163818c3aa7eSMark       ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr);
163918c3aa7eSMark     }
164018c3aa7eSMark   }
1641b7cbab4eSMark Adams   /* set options for subtype */
1642e55864a3SBarry Smith   if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
164318c3aa7eSMark 
164414a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
164514a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1646676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
16475b89ad90SMark F. Adams   PetscFunctionReturn(0);
16485b89ad90SMark F. Adams }
16495b89ad90SMark F. Adams 
16505b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
16515b89ad90SMark F. Adams /*MC
16521cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
16535b89ad90SMark F. Adams 
1654280d9858SJed Brown    Options Database Keys:
1655cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1656cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1657cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1658cab9ed1eSBarry 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
1659cab9ed1eSBarry 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>
1660cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1661cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
16626008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1663c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1664cab9ed1eSBarry Smith 
1665cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1666cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1667cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1668cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1669cab9ed1eSBarry Smith 
1670db9745e2SBarry Smith    Multigrid options:
1671db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1672db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1673db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1674cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
16755b89ad90SMark F. Adams 
16761cc46a46SBarry Smith 
167795452b02SPatrick Sanan   Notes:
167895452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1679db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1680db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1681db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
16821cc46a46SBarry Smith 
16835b89ad90SMark F. Adams   Level: intermediate
1684280d9858SJed Brown 
16851cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
168618c3aa7eSMark            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType()
16875b89ad90SMark F. Adams M*/
1688b2573a8aSBarry Smith 
16898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
16905b89ad90SMark F. Adams {
1691c1eae691SMark Adams   PetscErrorCode ierr,i;
16925b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
16935b89ad90SMark F. Adams   PC_MG          *mg;
16945b89ad90SMark F. Adams 
16955b89ad90SMark F. Adams   PetscFunctionBegin;
16961c1aac46SBarry Smith    /* register AMG type */
16971c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
16981c1aac46SBarry Smith 
16995b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
17001c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
17015b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
17025b89ad90SMark F. Adams 
17035b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1704b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
170569aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
17065b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
17075b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
17085b89ad90SMark F. Adams 
1709b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
17101ab5ffc9SJed Brown 
17119d5b6da9SMark F. Adams   pc_gamg->setup_count = 0;
17129d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
17139d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
17140a545947SLisandro Dalcin   pc_gamg->data    = NULL;
17155b89ad90SMark F. Adams 
17169d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
17175b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
17185b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
17195b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
17205b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
17215adeb434SBarry Smith   mg->view                = PCView_GAMG;
17225b89ad90SMark F. Adams 
172397d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
172497d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1725bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1726bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1727cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
172818c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr);
172918c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr);
173018c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr);
173118c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr);
17321cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1733cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1734171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1735ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1736ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1737bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1738c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1739bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1740c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1741bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
17429d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1743d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
17440c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1745171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1746a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1747a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1748038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
174925a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
175018c3aa7eSMark   for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1751c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
175218c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
17539ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1754d24ecf33SMark   ierr = PetscStrcpy(pc_gamg->esteig_type,NULL);CHKERRQ(ierr);
175518c3aa7eSMark   pc_gamg->esteig_max_it    = 10;
175618c3aa7eSMark   pc_gamg->use_sa_esteig    = -1;
175718c3aa7eSMark   pc_gamg->emin             = 0;
175818c3aa7eSMark   pc_gamg->emax             = 0;
175918c3aa7eSMark 
1760c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
17619d5b6da9SMark F. Adams 
1762bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1763bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
17645b89ad90SMark F. Adams   PetscFunctionReturn(0);
17655b89ad90SMark F. Adams }
17663e3471ccSMark Adams 
17673e3471ccSMark Adams /*@C
17683e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
17698a690491SBarry Smith     from PCInitializePackage().
17703e3471ccSMark Adams 
17713e3471ccSMark Adams  Level: developer
17723e3471ccSMark Adams 
17733e3471ccSMark Adams  .seealso: PetscInitialize()
17743e3471ccSMark Adams @*/
17753e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
17763e3471ccSMark Adams {
17773e3471ccSMark Adams   PetscErrorCode ierr;
17783e3471ccSMark Adams 
17793e3471ccSMark Adams   PetscFunctionBegin;
17803e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
17813e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
17823e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
17833e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
17848e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
17853e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1786c1c463dbSMark Adams 
1787c1c463dbSMark Adams   /* general events */
1788fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1789fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1790fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1791fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1792c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1793c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1794fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1795c1c463dbSMark Adams 
17965b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
17975b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
17985b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
17995b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
18005b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
18015b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
18025b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
18035b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
18045b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1805bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
18065b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
18075b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
18085b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
18095b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
18105b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
18115b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
18125b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
18135b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
18145b89ad90SMark F. Adams 
18155b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
18165b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
18175b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
18185b89ad90SMark F. Adams   /* create timer stages */
18195b89ad90SMark F. Adams #if defined GAMG_STAGES
18205b89ad90SMark F. Adams   {
18215b89ad90SMark F. Adams     char     str[32];
18225b89ad90SMark F. Adams     PetscInt lidx;
18235b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
18245b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
18255b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
18265b89ad90SMark F. Adams       sprintf(str,"MG Level %d",lidx);
18275b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
18285b89ad90SMark F. Adams     }
18295b89ad90SMark F. Adams   }
18305b89ad90SMark F. Adams #endif
18315b89ad90SMark F. Adams #endif
18323e3471ccSMark Adams   PetscFunctionReturn(0);
18333e3471ccSMark Adams }
18343e3471ccSMark Adams 
18353e3471ccSMark Adams /*@C
18361c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
18371c1aac46SBarry Smith     called from PetscFinalize() automatically.
18383e3471ccSMark Adams 
18393e3471ccSMark Adams  Level: developer
18403e3471ccSMark Adams 
18413e3471ccSMark Adams  .seealso: PetscFinalize()
18423e3471ccSMark Adams @*/
18433e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
18443e3471ccSMark Adams {
18453e3471ccSMark Adams   PetscErrorCode ierr;
18463e3471ccSMark Adams 
18473e3471ccSMark Adams   PetscFunctionBegin;
18483e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
18493e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
18503e3471ccSMark Adams   PetscFunctionReturn(0);
18513e3471ccSMark Adams }
1852a36cf38bSToby Isaac 
1853a36cf38bSToby Isaac /*@C
1854a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1855a36cf38bSToby Isaac 
1856a36cf38bSToby Isaac  Input Parameters:
1857a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1858a36cf38bSToby Isaac  - create - function for creating the gamg context.
1859a36cf38bSToby Isaac 
1860a36cf38bSToby Isaac   Level: advanced
1861a36cf38bSToby Isaac 
18621c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1863a36cf38bSToby Isaac @*/
1864a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1865a36cf38bSToby Isaac {
1866a36cf38bSToby Isaac   PetscErrorCode ierr;
1867a36cf38bSToby Isaac 
1868a36cf38bSToby Isaac   PetscFunctionBegin;
1869a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1870a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1871a36cf38bSToby Isaac   PetscFunctionReturn(0);
1872a36cf38bSToby Isaac }
1873a36cf38bSToby Isaac 
1874