xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision b4da3a1b2c3c2e292c3168fd4e6af360c056ffcc)
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 
9c9567895SMark #if defined(PETSC_HAVE_CUDA)
10c9567895SMark   #include <cuda_runtime.h>
11c9567895SMark #endif
12c9567895SMark 
13c9567895SMark #if defined(PETSC_HAVE_HIP)
14c9567895SMark   #include <hip/hip_runtime.h>
15c9567895SMark #endif
16c9567895SMark 
170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
180cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
19b4fbaa2aSMark F. Adams #endif
200cbbd2e1SMark F. Adams 
210cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
22fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
23fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
240cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
250cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
260cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
270cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
28fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
290cbbd2e1SMark F. Adams #endif
300cbbd2e1SMark F. Adams 
31b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
320cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
3318c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
34b4fbaa2aSMark F. Adams #endif
35f96513f1SMatthew G Knepley 
360a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
373e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
389d5b6da9SMark F. Adams 
39d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
40d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
41d3d6bff4SMark F. Adams {
4218c3aa7eSMark   PetscErrorCode ierr, level;
43d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
44d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
45d3d6bff4SMark F. Adams 
46d3d6bff4SMark F. Adams   PetscFunctionBegin;
4722a233eaSStefano Zampini   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
481c1aac46SBarry Smith   pc_gamg->data_sz = 0;
49878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
5018c3aa7eSMark   for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
5118c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
5218c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
5318c3aa7eSMark   }
5418c3aa7eSMark   pc_gamg->emin = 0;
5518c3aa7eSMark   pc_gamg->emax = 0;
56a2f3521dSMark F. Adams   PetscFunctionReturn(0);
57a2f3521dSMark F. Adams }
58a2f3521dSMark F. Adams 
595b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
605b89ad90SMark F. Adams /*
61c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
62a147abb0SMark F. Adams      of active processors.
635b89ad90SMark F. Adams 
645b89ad90SMark F. Adams    Input Parameter:
65a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
66a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
679d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
68c5bfad50SMark F. Adams    . cr_bs - coarse block size
693530afc2SMark F. Adams    In/Output Parameter:
70a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
71afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
7211e60469SMark F. Adams    Output Parameter:
733530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
745b89ad90SMark F. Adams */
755cb416c2SMark F. Adams 
76171cca9aSMark 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)
775b89ad90SMark F. Adams {
78a2f3521dSMark F. Adams   PetscErrorCode  ierr;
799d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
80486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
81a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
823b4367a7SBarry Smith   MPI_Comm        comm;
83c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
843ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
855b89ad90SMark F. Adams 
865b89ad90SMark F. Adams   PetscFunctionBegin;
873b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
88ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
89ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
90c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
919d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
92038e3b61SMark F. Adams 
93ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
94ce7c7f2fSMark Adams 
953ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
960298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
973ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
983ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9973911c69SBarry Smith   } else {
1003ae0bb68SMark Adams     PetscInt  bs;
1013ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
1023ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
1033ae0bb68SMark Adams   }
104c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
105c9567895SMark   if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level==0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */
106c9567895SMark #if defined(PETSC_HAVE_CUDA)
107c9567895SMark     PetscShmComm pshmcomm;
108c9567895SMark     PetscMPIInt  locrank;
109c9567895SMark     MPI_Comm     loccomm;
110c9567895SMark     PetscInt     s_nnodes,r_nnodes, new_new_size;
111c9567895SMark     cudaError_t  cerr;
112c9567895SMark     int          devCount;
113c9567895SMark     ierr = PetscShmCommGet(comm,&pshmcomm);CHKERRQ(ierr);
114c9567895SMark     ierr = PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);CHKERRQ(ierr);
115c9567895SMark     ierr = MPI_Comm_rank(loccomm, &locrank);CHKERRQ(ierr);
116c9567895SMark     s_nnodes = !locrank;
117c9567895SMark     ierr = MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
118c9567895SMark     if (size%r_nnodes) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes);
119c9567895SMark     devCount = 0;
120c9567895SMark     cerr = cudaGetDeviceCount(&devCount);
121c9567895SMark     cudaGetLastError(); /* Reset the last error */
122c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
123c9567895SMark       new_new_size = r_nnodes * devCount;
124c9567895SMark       new_size = new_new_size;
125c9567895SMark       ierr = PetscInfo5(pc,"Fine grid with Cuda. %D nodes. Change new active set size %D --> %D (devCount=%D #nodes=%D)\n",r_nnodes,nactive,new_size,devCount,r_nnodes);CHKERRQ(ierr);
126c9567895SMark     } else {
127c9567895SMark       ierr = PetscInfo(pc,"With Cuda but no device. Use heuristics.");CHKERRQ(ierr);
128c9567895SMark       goto HEURISTIC;
129c9567895SMark     }
130c9567895SMark #else
131c9567895SMark     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
132c9567895SMark #endif
133c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
134c9567895SMark     if (nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of active process %D wrt reduction factor %D",nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
135c9567895SMark     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
136c9567895SMark     ierr = PetscInfo3(pc,"Manually setting reduction to %D active processes (%D/%D)\n",new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);CHKERRQ(ierr);
137c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
138c9567895SMark     new_size = 1;
1398abdc6daSStefano Zampini     ierr = PetscInfo1(pc,"Force coarsest grid reduction to %D active processes\n",new_size);CHKERRQ(ierr);
140c9567895SMark   } else {
141472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
142c9567895SMark #if defined(PETSC_HAVE_CUDA)
143c9567895SMark     HEURISTIC:
144c9567895SMark #endif
1450298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
146a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1477f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
148c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
1498abdc6daSStefano Zampini     ierr = PetscInfo2(pc,"Coarse grid reduction from %D to %D active processes\n",nactive,new_size);CHKERRQ(ierr);
150a2f3521dSMark F. Adams   }
1512e3501ffSMark Adams   if (new_size==nactive) {
152ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
153ce7c7f2fSMark Adams     if (new_size < size) {
154ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
155c9567895SMark       ierr = PetscInfo1(pc,"reduced grid using same number of processors (%D) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr);
156ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
157b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
158b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
159ce7c7f2fSMark Adams       }
160ce7c7f2fSMark Adams     }
161ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1622e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
163192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
164885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
16571959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
16671959b99SBarry 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);
167ce7c7f2fSMark Adams     /* get new_size and rfactor */
168ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
169ce7c7f2fSMark Adams       /* find factor */
170ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
171ce7c7f2fSMark Adams       else {
172ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
173ce7c7f2fSMark Adams         jj = -1;
174ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
175ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
176ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
177ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
178ce7c7f2fSMark Adams             if (fact > best_fact) {
179ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
180ce7c7f2fSMark Adams             }
181ce7c7f2fSMark Adams           }
182ce7c7f2fSMark Adams         }
183ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
184ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
185ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
186ce7c7f2fSMark Adams         else expand_factor = rfactor;
187ce7c7f2fSMark Adams       }
188ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
1894cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1904cdfd227SMark         *a_Amat_crs = Cmat;
191c9567895SMark         ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr);
192ce7c7f2fSMark Adams         PetscFunctionReturn(0);
193ce7c7f2fSMark Adams       }
194ce7c7f2fSMark Adams     }
1954cdfd227SMark #if defined PETSC_GAMG_USE_LOG
1964cdfd227SMark     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
1974cdfd227SMark #endif
198a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
199785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
2002e3501ffSMark Adams     if (pc_gamg->repart) {
201a2f3521dSMark F. Adams       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
2025a9b9e01SMark F. Adams       Mat      adj;
203c9567895SMark       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);
204a2f3521dSMark F. Adams       /* get 'adj' */
205c5bfad50SMark F. Adams       if (cr_bs == 1) {
206038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
207806fa848SBarry Smith       } else {
208a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
209eb07cef2SMark F. Adams         Mat               tMat;
210a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
211b4fbaa2aSMark F. Adams         const PetscScalar *vals;
212b4fbaa2aSMark F. Adams         const PetscInt    *idx;
213a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
21439d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
215d9558ea9SBarry Smith         MatType           mtype;
216b4fbaa2aSMark F. Adams 
217e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
218a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
219a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
220c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
2210a545947SLisandro Dalcin           ierr      = MatGetRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
222c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
223c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
2240a545947SLisandro Dalcin           ierr      = MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
2253ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
2263ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
22758471d46SMark F. Adams         }
2286876a03eSMark F. Adams 
229d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
2303b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
2313ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
232d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
233a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
234a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
235e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
236eb07cef2SMark F. Adams 
237a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
238c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
23922063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
240eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
241c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
242eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
243eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
244eb07cef2SMark F. Adams           }
24522063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
246eb07cef2SMark F. Adams         }
247eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249eb07cef2SMark F. Adams 
250b4fbaa2aSMark F. Adams         if (llev++ == -1) {
251b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2528caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2533b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
254b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2553bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
256b4fbaa2aSMark F. Adams         }
257eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
258eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
259a2f3521dSMark F. Adams       } /* create 'adj' */
260f150b916SMark F. Adams 
261a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2625a9b9e01SMark F. Adams         char            prefix[256];
2635a9b9e01SMark F. Adams         const char      *pcpre;
264b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
265b4fbaa2aSMark F. Adams         MatPartitioning mpart;
266a4b7d37bSMark F. Adams         IS              proc_is;
2672f03bc48SMark F. Adams 
2683b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2695ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2709d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2718caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
27259a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
27311e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
274c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
275a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27611e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2775a9b9e01SMark F. Adams 
2785ef31b24SMark F. Adams         /* collect IS info */
279785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
280a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
281a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
282c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
283ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
284eb07cef2SMark F. Adams           }
2855ef31b24SMark F. Adams         }
286a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
287a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2885ef31b24SMark F. Adams       }
2895ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2905a9b9e01SMark F. Adams 
2913b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2928263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
29331cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
294ce7c7f2fSMark Adams       PetscInt targetPE;
2954cdfd227SMark       if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
296302440fdSBarry Smith       ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr);
297ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2983b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
299a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
300e33ef3b1SMark F. Adams 
30111e60469SMark F. Adams     /*
302a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
30311e60469SMark F. Adams     */
304a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
3057700e67bSMark Adams     is_eq_num_prim = is_eq_num;
30611e60469SMark F. Adams     /*
307a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
30811e60469SMark F. Adams     */
309c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
310c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
311a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
312ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
313a2f3521dSMark F. Adams 
314a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
315885364a3SMark 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 */
316885364a3SMark Adams     {
317885364a3SMark Adams       Vec            src_crd, dest_crd;
318885364a3SMark 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;
319885364a3SMark Adams       VecScatter     vecscat;
320885364a3SMark Adams       PetscScalar    *array;
321885364a3SMark Adams       IS isscat;
322a2f3521dSMark F. Adams       /* move data (for primal equations only) */
32322063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
3243b4367a7SBarry Smith       ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
3253ae0bb68SMark Adams       ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
326c0dedaeaSBarry Smith       ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
32711e60469SMark F. Adams       /*
3289d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
329c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
33011e60469SMark F. Adams       */
331854ce69bSBarry Smith       ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
332a2f3521dSMark F. Adams       ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3333ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
334c5bfad50SMark F. Adams         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
335a2f3521dSMark F. Adams         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
33611e60469SMark F. Adams       }
337a2f3521dSMark F. Adams       ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3383ae0bb68SMark Adams       ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
33992a756f0SMark F. Adams       ierr = PetscFree(tidx);CHKERRQ(ierr);
34011e60469SMark F. Adams       /*
34111e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
34211e60469SMark F. Adams       */
3433ae0bb68SMark Adams       ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3449d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3453ae0bb68SMark Adams         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3463ae0bb68SMark Adams         for (ii=0; ii<ncrs; ii++) {
3479d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
348a2f3521dSMark F. Adams             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
349c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
350676e1743SMark F. Adams             ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
351d3d6bff4SMark F. Adams           }
352038e3b61SMark F. Adams         }
353eb07cef2SMark F. Adams       }
354eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
355eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
35611e60469SMark F. Adams       /*
35711e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
35811e60469SMark F. Adams         to the correct processor
35911e60469SMark F. Adams       */
3609448b7f1SJunchao Zhang       ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
36111e60469SMark F. Adams       ierr = ISDestroy(&isscat);CHKERRQ(ierr);
36211e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
36311e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
36411e60469SMark F. Adams       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
36511e60469SMark F. Adams       ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
36611e60469SMark F. Adams       /*
36711e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
36811e60469SMark F. Adams       */
369c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
370578f55a3SPeter Brune       ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3712fa5cd67SKarl Rupp 
3723ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3733ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3742fa5cd67SKarl Rupp 
37511e60469SMark F. Adams       ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3769d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3773ae0bb68SMark Adams         for (ii=0; ii<ncrs_new; ii++) {
3789d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
379a2f3521dSMark F. Adams             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
380c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
381d3d6bff4SMark F. Adams           }
382038e3b61SMark F. Adams         }
383038e3b61SMark F. Adams       }
38411e60469SMark F. Adams       ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
38511e60469SMark F. Adams       ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
386885364a3SMark Adams     }
387a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3880cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3890cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
390ed3f9983SMark F. Adams #endif
39111e60469SMark F. Adams     /*
3927dae84e0SHong Zhang       Invert for MatCreateSubMatrix
39311e60469SMark F. Adams     */
394a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
395a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
396c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
397a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
398a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
399a2f3521dSMark F. Adams     }
4003cb8563fSToby Isaac     if (Pcolumnperm) {
4013cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
4023cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
4033cb8563fSToby Isaac     }
404a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
4050cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4060cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
4070cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
408ed3f9983SMark F. Adams #endif
409a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
410a2f3521dSMark F. Adams     {
411a2f3521dSMark F. Adams       Mat mat;
4127dae84e0SHong Zhang       ierr        = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
413a2f3521dSMark F. Adams       *a_Amat_crs = mat;
414a2f3521dSMark F. Adams     }
415038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
416a2f3521dSMark F. Adams 
4170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4180cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
419ed3f9983SMark F. Adams #endif
42011e60469SMark F. Adams     /* prolongator */
42111e60469SMark F. Adams     {
42211e60469SMark F. Adams       IS       findices;
423a2f3521dSMark F. Adams       PetscInt Istart,Iend;
424a2f3521dSMark F. Adams       Mat      Pnew;
42562294041SBarry Smith 
426a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4270cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4280cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
429ed3f9983SMark F. Adams #endif
4303b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
431c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
4327dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
43311e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
434e589036eSStefano Zampini #if defined(PETSC_HAVE_CUDA)
435e589036eSStefano Zampini       ierr = MatAIJCUSPARSESetGenerateTranspose(Pnew,PETSC_TRUE);CHKERRQ(ierr);
436e589036eSStefano Zampini #endif
437c5bfad50SMark F. Adams 
4380cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4390cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
440ed3f9983SMark F. Adams #endif
4413530afc2SMark F. Adams       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
442a2f3521dSMark F. Adams 
443a2f3521dSMark F. Adams       /* output - repartitioned */
444a2f3521dSMark F. Adams       *a_P_inout = Pnew;
445e33ef3b1SMark F. Adams     }
446a2f3521dSMark F. Adams     ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr);
4475b89ad90SMark F. Adams 
448c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
449ce7c7f2fSMark Adams 
450ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
451ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
452ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4538bca76a6SMark Adams       static PetscInt llev = 2;
45439d09545SMark Adams       ierr = PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr);
455ce7c7f2fSMark Adams #endif
456b470e4b4SRichard Tran Mills       ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
457b470e4b4SRichard Tran Mills       ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
458ce7c7f2fSMark Adams       if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
459ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
460ce7c7f2fSMark Adams         PetscMPIInt size;
461ffc4695bSBarry Smith         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
462ce7c7f2fSMark Adams         if (size > 1) {
463ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
464b470e4b4SRichard Tran Mills           ierr = VecBindToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr);
465b470e4b4SRichard Tran Mills           ierr = VecBindToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr);
466ce7c7f2fSMark Adams         }
467ce7c7f2fSMark Adams       }
468ce7c7f2fSMark Adams     }
4694cdfd227SMark #if defined PETSC_GAMG_USE_LOG
4704cdfd227SMark     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
4714cdfd227SMark #endif
472a2f3521dSMark F. Adams   }
4735b89ad90SMark F. Adams   PetscFunctionReturn(0);
4745b89ad90SMark F. Adams }
4755b89ad90SMark F. Adams 
4764b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
4774b1575e2SStefano Zampini {
4784b1575e2SStefano Zampini   PetscErrorCode ierr;
4794b1575e2SStefano Zampini   const char     *prefix;
4804b1575e2SStefano Zampini   char           addp[32];
4814b1575e2SStefano Zampini   PC_MG          *mg      = (PC_MG*)a_pc->data;
4824b1575e2SStefano Zampini   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4834b1575e2SStefano Zampini 
4844b1575e2SStefano Zampini   PetscFunctionBegin;
4854b1575e2SStefano Zampini   ierr = PCGetOptionsPrefix(a_pc,&prefix);CHKERRQ(ierr);
4864b1575e2SStefano Zampini   ierr = PetscInfo1(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);CHKERRQ(ierr);
4874b1575e2SStefano Zampini   ierr = MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);CHKERRQ(ierr);
4884b1575e2SStefano Zampini   ierr = MatSetOptionsPrefix(*Gmat2,prefix);CHKERRQ(ierr);
4894b1575e2SStefano Zampini   ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);CHKERRQ(ierr);
4904b1575e2SStefano Zampini   ierr = MatAppendOptionsPrefix(*Gmat2,addp);CHKERRQ(ierr);
491*b4da3a1bSStefano Zampini   if ((*Gmat2)->structurally_symmetric) {
492*b4da3a1bSStefano Zampini     ierr = MatProductSetType(*Gmat2,MATPRODUCT_AB);CHKERRQ(ierr);
493*b4da3a1bSStefano Zampini   } else {
494*b4da3a1bSStefano Zampini #if defined(PETSC_HAVE_CUDA)
495*b4da3a1bSStefano Zampini     ierr = MatAIJCUSPARSESetGenerateTranspose(Gmat1,PETSC_TRUE);CHKERRQ(ierr);
496*b4da3a1bSStefano Zampini #endif
4974b1575e2SStefano Zampini     ierr = MatProductSetType(*Gmat2,MATPRODUCT_AtB);CHKERRQ(ierr);
498*b4da3a1bSStefano Zampini   }
4994b1575e2SStefano Zampini   ierr = MatProductSetFromOptions(*Gmat2);CHKERRQ(ierr);
5004b1575e2SStefano Zampini   ierr = MatProductSymbolic(*Gmat2);CHKERRQ(ierr);
501*b4da3a1bSStefano Zampini   ierr = MatProductClear(*Gmat2);CHKERRQ(ierr);
5024b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
5034b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
5044b1575e2SStefano Zampini   PetscFunctionReturn(0);
5054b1575e2SStefano Zampini }
5064b1575e2SStefano Zampini 
5075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5085b89ad90SMark F. Adams /*
5095b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5105b89ad90SMark F. Adams                     by setting data structures and options.
5115b89ad90SMark F. Adams 
5125b89ad90SMark F. Adams    Input Parameter:
5135b89ad90SMark F. Adams .  pc - the preconditioner context
5145b89ad90SMark F. Adams 
5155b89ad90SMark F. Adams */
5169d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5175b89ad90SMark F. Adams {
5185b89ad90SMark F. Adams   PetscErrorCode ierr;
5199d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5205b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5212adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
52218c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
5233b4367a7SBarry Smith   MPI_Comm       comm;
524c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
52518c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
52618c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
527a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
528569f4572SMark Adams   MatInfo        info;
529171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
5305ef31b24SMark F. Adams 
5315b89ad90SMark F. Adams   PetscFunctionBegin;
5323b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
533ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
534ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
535dfd5c07aSMark F. Adams 
5368abdc6daSStefano Zampini   if (pc->setupcalled) {
5378abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
538878e152fSMark F. Adams       /* reset everything */
539878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
540878e152fSMark F. Adams       pc->setupcalled = 0;
541806fa848SBarry Smith     } else {
54284d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
54303a628feSMark F. Adams       /* just do Galerkin grids */
54458471d46SMark F. Adams       Mat          B,dA,dB;
54558471d46SMark F. Adams 
5469d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
54758471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
54823ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
54958471d46SMark F. Adams         /* (re)set to get dirty flag */
55023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
55158471d46SMark F. Adams 
5522fb0b348SMark F. Adams         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
5538abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX ;
5548abdc6daSStefano Zampini 
5558abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
55623ee1639SBarry Smith           ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
5578abdc6daSStefano Zampini           if (B->product) {
5588abdc6daSStefano Zampini             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
5598abdc6daSStefano Zampini               reuse = MAT_REUSE_MATRIX;
56003a628feSMark F. Adams             }
5618abdc6daSStefano Zampini           }
5628abdc6daSStefano Zampini           if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); }
5638abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
5648abdc6daSStefano Zampini             ierr = PetscInfo1(pc,"RAP after first solve, reuse matrix level %D\n",level);CHKERRQ(ierr);
5658abdc6daSStefano Zampini           } else {
5668abdc6daSStefano Zampini             ierr = PetscInfo1(pc,"RAP after first solve, new matrix level %D\n",level);CHKERRQ(ierr);
5678abdc6daSStefano Zampini           }
5688abdc6daSStefano Zampini           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);CHKERRQ(ierr);
5698abdc6daSStefano Zampini           mglevels[level]->A = B;
57023ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
57158471d46SMark F. Adams           dB   = B;
57258471d46SMark F. Adams         }
5735f8cf99dSMark F. Adams       }
574d5280255SMark F. Adams 
575d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
57658471d46SMark F. Adams       PetscFunctionReturn(0);
577eb07cef2SMark F. Adams     }
578878e152fSMark F. Adams   }
579f6536408SMark F. Adams 
580878e152fSMark F. Adams   if (!pc_gamg->data) {
581878e152fSMark F. Adams     if (pc_gamg->orig_data) {
582878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5830298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5842fa5cd67SKarl Rupp 
585878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
586878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
587878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5882fa5cd67SKarl Rupp 
589785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
590878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
591806fa848SBarry Smith     } else {
5921ab5ffc9SJed Brown       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5937700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5949d5b6da9SMark F. Adams     }
595878e152fSMark F. Adams   }
596878e152fSMark F. Adams 
597878e152fSMark F. Adams   /* cache original data for reuse */
5981c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
599785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
600878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
601878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
602878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
603878e152fSMark F. Adams   }
604038e3b61SMark F. Adams 
605302f38e8SMark F. Adams   /* get basic dims */
606302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
607171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
60884d3f75bSMark F. Adams 
609569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
610569f4572SMark Adams   nnz0   = info.nz_used;
611569f4572SMark Adams   nnztot = info.nz_used;
612c9567895SMark   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);
613569f4572SMark Adams 
614a2f3521dSMark F. Adams   /* Get A_i and R_i */
61562294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
6169ab59c8bSMark Adams     pc_gamg->current_level = level;
61718c3aa7eSMark     if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
6185b89ad90SMark F. Adams     level1 = level + 1;
6190cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6200cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
621a2f3521dSMark F. Adams #if (defined GAMG_STAGES)
622a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
623b4fbaa2aSMark F. Adams #endif
624a2f3521dSMark F. Adams #endif
625c8b0795cSMark F. Adams     { /* construct prolongator */
626725b86d8SJed Brown       Mat              Gmat;
6270cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6287700e67bSMark Adams       Mat              Prol11;
629c8b0795cSMark F. Adams 
6307700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
6311ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
6327700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
633c8b0795cSMark F. Adams 
634a2f3521dSMark F. Adams       /* could have failed to create new level */
635a2f3521dSMark F. Adams       if (Prol11) {
636f7df55f0SStefano Zampini         const char *prefix;
637f7df55f0SStefano Zampini         char       addp[32];
638f7df55f0SStefano Zampini 
6399d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6400298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
641a2f3521dSMark F. Adams 
642fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
643c8b0795cSMark F. Adams           /* smooth */
644fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
645c8b0795cSMark F. Adams         }
646c8b0795cSMark F. Adams 
6470c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
6481b18a24aSMark Adams           PetscInt bs;
6491b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
6500a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
651ffc955d6SMark F. Adams         }
652ffc955d6SMark F. Adams 
653f7df55f0SStefano Zampini         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
654f7df55f0SStefano Zampini         ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr);
655c9567895SMark         ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);CHKERRQ(ierr);
656f7df55f0SStefano Zampini         ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr);
65791f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
658f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
659f7df55f0SStefano Zampini #if defined(PETSC_HAVE_CUDA)
660e589036eSStefano Zampini         ierr = MatAIJCUSPARSESetGenerateTranspose(Prol11,PETSC_TRUE);CHKERRQ(ierr);
661f7df55f0SStefano Zampini #endif
662f7df55f0SStefano Zampini         ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr);
6634bde40a0SMark Adams         Parr[level1] = Prol11;
6644bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6654bde40a0SMark Adams 
666a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
66741b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
668a2f3521dSMark F. Adams     } /* construct prolongator scope */
6690cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6700cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
671c8b0795cSMark F. Adams #endif
6727f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
673171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
674569f4572SMark Adams       ierr =  PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr);
67562294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
676a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
677a90e85d9SMark Adams #endif
678c8b0795cSMark F. Adams       break;
679c8b0795cSMark F. Adams     }
6800cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6810cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
682b4fbaa2aSMark F. Adams #endif
683171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
684c9567895SMark     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
685171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6860e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
687171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
688a2f3521dSMark F. Adams 
6890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
6900cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
691b4fbaa2aSMark F. Adams #endif
692171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
693569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
694569f4572SMark Adams     nnztot += info.nz_used;
695c9567895SMark     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);
696569f4572SMark Adams 
6970cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
698b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
699b4fbaa2aSMark F. Adams #endif
700a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
7019ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
7029ab59c8bSMark Adams       ierr =  PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr);
703a90e85d9SMark Adams       level++;
704a90e85d9SMark Adams       break;
705a90e85d9SMark Adams     }
706c8b0795cSMark F. Adams   } /* levels */
707c8b0795cSMark F. Adams   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
708c8b0795cSMark F. Adams 
709569f4572SMark Adams   ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr);
7109d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7115b89ad90SMark F. Adams   fine_level       = level;
7120298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
7135b89ad90SMark F. Adams 
71462294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7150ed2132dSStefano Zampini     PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP);
7160ed2132dSStefano Zampini 
717d5280255SMark F. Adams     /* set default smoothers & set operators */
71862294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
719ffc955d6SMark F. Adams       KSP smoother;
720ffc955d6SMark F. Adams       PC  subpc;
721a2f3521dSMark F. Adams 
7229d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
723f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
724ffc955d6SMark F. Adams 
725a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
726a2f3521dSMark F. Adams       /* set ops */
72723ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
728a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
729a2f3521dSMark F. Adams 
730a2f3521dSMark F. Adams       /* set defaults */
7316c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
732a2f3521dSMark F. Adams 
7330c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
7340c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
7352d3561bbSSatish Balay         PetscInt sz;
7367a28f3e5SMark Adams         IS       *iss;
737a2f3521dSMark F. Adams 
7382d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7397a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
7400c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
7410a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
7420c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
7437f66b68fSMark Adams         if (!sz) {
744ffc955d6SMark F. Adams           IS       is;
7450a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7467a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
747a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
748806fa848SBarry Smith         } else {
749a94c3b12SMark F. Adams           PetscInt kk;
7507a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
751a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
7527a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
753a94c3b12SMark F. Adams           }
7547a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
755ffc955d6SMark F. Adams         }
7560298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
757ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
758806fa848SBarry Smith       } else {
759890ffe84SMark Adams         ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr);
760ffc955d6SMark F. Adams       }
761d5280255SMark F. Adams     }
762d5280255SMark F. Adams     {
763d5280255SMark F. Adams       /* coarse grid */
764d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
765d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
7660ed2132dSStefano Zampini 
767d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
76823ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
769cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
770d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
771d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
772d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
773d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
77471959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
77571959b99SBarry Smith         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
776d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
777d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7789dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7792fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
78008e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
7815b42dca8SJed Brown         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
7825b42dca8SJed Brown          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
7835b42dca8SJed Brown          * view every subdomain as though they were different. */
7845b42dca8SJed Brown         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
785d5280255SMark F. Adams       }
786cf8ae1d3SMark Adams     }
787d5280255SMark F. Adams 
788d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
789d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
790e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
791d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
79269aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
793d5280255SMark F. Adams 
79418c3aa7eSMark     /* setup cheby eigen estimates from SA */
7950ed2132dSStefano Zampini     if (pc_gamg->use_sa_esteig==1) {
79618c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
79718c3aa7eSMark         KSP       smoother;
79818c3aa7eSMark         PetscBool ischeb;
7990ed2132dSStefano Zampini 
8000ed2132dSStefano Zampini         savesetfromoptions[level] = NULL;
80118c3aa7eSMark         ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
80218c3aa7eSMark         ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr);
80318c3aa7eSMark         if (ischeb) {
80418c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
8050ed2132dSStefano Zampini 
8060ed2132dSStefano Zampini           ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr); /* let command line emax override using SA's eigenvalues */
8070ed2132dSStefano Zampini           if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) {
80818c3aa7eSMark             PC        subpc;
80918c3aa7eSMark             PetscBool isjac;
81018c3aa7eSMark             ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
81118c3aa7eSMark             ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr);
8120ed2132dSStefano Zampini             if (isjac && pc_gamg->use_sa_esteig==1) {
81318c3aa7eSMark               PetscReal emax,emin;
8140ed2132dSStefano Zampini 
81518c3aa7eSMark               emin = mg->min_eigen_DinvA[level];
81618c3aa7eSMark               emax = mg->max_eigen_DinvA[level];
81718c3aa7eSMark               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);
81818c3aa7eSMark               cheb->emin_computed = emin;
81918c3aa7eSMark               cheb->emax_computed = emax;
82018c3aa7eSMark               ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr);
8210ed2132dSStefano Zampini 
8220ed2132dSStefano Zampini               /* We have set the eigenvalues and consumed the transformation values
8230ed2132dSStefano Zampini                  prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG
8240ed2132dSStefano Zampini                  below when setfromoptions will be called again */
8250ed2132dSStefano Zampini               savesetfromoptions[level] = smoother->ops->setfromoptions;
8260ed2132dSStefano Zampini               smoother->ops->setfromoptions = NULL;
82718c3aa7eSMark             }
82818c3aa7eSMark           }
82918c3aa7eSMark         }
83018c3aa7eSMark       }
8310ed2132dSStefano Zampini     }
8320ed2132dSStefano Zampini 
8330ed2132dSStefano Zampini     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8340ed2132dSStefano Zampini 
8350ed2132dSStefano Zampini     /* restore Chebyshev smoother for next calls */
8360ed2132dSStefano Zampini     if (pc_gamg->use_sa_esteig==1) {
8370ed2132dSStefano Zampini       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
8380ed2132dSStefano Zampini         if (savesetfromoptions[level]) {
8390ed2132dSStefano Zampini           KSP smoother;
8400ed2132dSStefano Zampini           ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
8410ed2132dSStefano Zampini           smoother->ops->setfromoptions = savesetfromoptions[level];
8420ed2132dSStefano Zampini         }
8430ed2132dSStefano Zampini       }
8440ed2132dSStefano Zampini     }
84518c3aa7eSMark 
846d5280255SMark F. Adams     /* clean up */
847d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
848587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
849587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8505b89ad90SMark F. Adams     }
851806fa848SBarry Smith   } else {
8525f8cf99dSMark F. Adams     KSP smoother;
8530ed2132dSStefano Zampini 
854302440fdSBarry Smith     ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
8559d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
85623ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8575f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8589d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8595f8cf99dSMark F. Adams   }
8605b89ad90SMark F. Adams   PetscFunctionReturn(0);
8615b89ad90SMark F. Adams }
8625b89ad90SMark F. Adams 
863eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8645b89ad90SMark F. Adams /*
8655b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8665b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8675b89ad90SMark F. Adams 
8685b89ad90SMark F. Adams    Input Parameter:
8695b89ad90SMark F. Adams .  pc - the preconditioner context
8705b89ad90SMark F. Adams 
8715b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8725b89ad90SMark F. Adams */
8735b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8745b89ad90SMark F. Adams {
8755b89ad90SMark F. Adams   PetscErrorCode ierr;
8765b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8775b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8785b89ad90SMark F. Adams 
8795b89ad90SMark F. Adams   PetscFunctionBegin;
8805b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8819b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8829b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
8839b8ffb57SJed Brown   }
8841ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
8851ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
8865b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8875b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8885b89ad90SMark F. Adams   PetscFunctionReturn(0);
8895b89ad90SMark F. Adams }
8905b89ad90SMark F. Adams 
891676e1743SMark F. Adams /*@
892cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
893676e1743SMark F. Adams 
8941cc46a46SBarry Smith    Logically Collective on PC
895676e1743SMark F. Adams 
896676e1743SMark F. Adams    Input Parameters:
8971cc46a46SBarry Smith +  pc - the preconditioner context
8981cc46a46SBarry Smith -  n - the number of equations
899676e1743SMark F. Adams 
900676e1743SMark F. Adams 
901676e1743SMark F. Adams    Options Database Key:
9021cc46a46SBarry Smith .  -pc_gamg_process_eq_limit <limit>
903676e1743SMark F. Adams 
90495452b02SPatrick Sanan    Notes:
90595452b02SPatrick 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
906cab9ed1eSBarry Smith           that has degrees of freedom
907cab9ed1eSBarry Smith 
908676e1743SMark F. Adams    Level: intermediate
909676e1743SMark F. Adams 
910c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
911676e1743SMark F. Adams @*/
912676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
913676e1743SMark F. Adams {
914676e1743SMark F. Adams   PetscErrorCode ierr;
915676e1743SMark F. Adams 
916676e1743SMark F. Adams   PetscFunctionBegin;
917676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
918676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
919676e1743SMark F. Adams   PetscFunctionReturn(0);
920676e1743SMark F. Adams }
921676e1743SMark F. Adams 
9221e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
923676e1743SMark F. Adams {
924c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
925c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
926676e1743SMark F. Adams 
927676e1743SMark F. Adams   PetscFunctionBegin;
9289d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
929676e1743SMark F. Adams   PetscFunctionReturn(0);
930676e1743SMark F. Adams }
931676e1743SMark F. Adams 
932389730f3SMark F. Adams /*@
933cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
934389730f3SMark F. Adams 
935389730f3SMark F. Adams  Collective on PC
936389730f3SMark F. Adams 
937389730f3SMark F. Adams    Input Parameters:
9381cc46a46SBarry Smith +  pc - the preconditioner context
9391cc46a46SBarry Smith -  n - maximum number of equations to aim for
940389730f3SMark F. Adams 
941389730f3SMark F. Adams    Options Database Key:
9421cc46a46SBarry Smith .  -pc_gamg_coarse_eq_limit <limit>
943389730f3SMark F. Adams 
94474329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
94574329af1SBarry Smith      has less than 1000 unknowns.
94674329af1SBarry Smith 
947389730f3SMark F. Adams    Level: intermediate
948389730f3SMark F. Adams 
949c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
950389730f3SMark F. Adams @*/
951389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
952389730f3SMark F. Adams {
953389730f3SMark F. Adams   PetscErrorCode ierr;
954389730f3SMark F. Adams 
955389730f3SMark F. Adams   PetscFunctionBegin;
956389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
957389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
958389730f3SMark F. Adams   PetscFunctionReturn(0);
959389730f3SMark F. Adams }
960389730f3SMark F. Adams 
9611e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
962389730f3SMark F. Adams {
963389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
964389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
965389730f3SMark F. Adams 
966389730f3SMark F. Adams   PetscFunctionBegin;
9679d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
968389730f3SMark F. Adams   PetscFunctionReturn(0);
969389730f3SMark F. Adams }
970389730f3SMark F. Adams 
971676e1743SMark F. Adams /*@
972cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
973676e1743SMark F. Adams 
974676e1743SMark F. Adams    Collective on PC
975676e1743SMark F. Adams 
976676e1743SMark F. Adams    Input Parameters:
9771cc46a46SBarry Smith +  pc - the preconditioner context
9781cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
979676e1743SMark F. Adams 
980676e1743SMark F. Adams    Options Database Key:
9811cc46a46SBarry Smith .  -pc_gamg_repartition <true,false>
982676e1743SMark F. Adams 
98395452b02SPatrick Sanan    Notes:
98495452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
985cab9ed1eSBarry Smith 
986676e1743SMark F. Adams    Level: intermediate
987676e1743SMark F. Adams 
988676e1743SMark F. Adams .seealso: ()
989676e1743SMark F. Adams @*/
990cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
991676e1743SMark F. Adams {
992676e1743SMark F. Adams   PetscErrorCode ierr;
993676e1743SMark F. Adams 
994676e1743SMark F. Adams   PetscFunctionBegin;
995676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
996cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
997676e1743SMark F. Adams   PetscFunctionReturn(0);
998676e1743SMark F. Adams }
999676e1743SMark F. Adams 
1000cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
1001676e1743SMark F. Adams {
1002c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1003c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1004676e1743SMark F. Adams 
1005676e1743SMark F. Adams   PetscFunctionBegin;
10069d5b6da9SMark F. Adams   pc_gamg->repart = n;
1007676e1743SMark F. Adams   PetscFunctionReturn(0);
1008676e1743SMark F. Adams }
1009676e1743SMark F. Adams 
1010dfd5c07aSMark F. Adams /*@
101118c3aa7eSMark    PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby)
101218c3aa7eSMark 
101318c3aa7eSMark    Collective on PC
101418c3aa7eSMark 
101518c3aa7eSMark    Input Parameters:
101618c3aa7eSMark +  pc - the preconditioner context
101718c3aa7eSMark -  n - number of its
101818c3aa7eSMark 
101918c3aa7eSMark    Options Database Key:
102018c3aa7eSMark .  -pc_gamg_esteig_ksp_max_it <its>
102118c3aa7eSMark 
102218c3aa7eSMark    Notes:
102318c3aa7eSMark 
102418c3aa7eSMark    Level: intermediate
102518c3aa7eSMark 
102618c3aa7eSMark .seealso: ()
102718c3aa7eSMark @*/
102818c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n)
102918c3aa7eSMark {
103018c3aa7eSMark   PetscErrorCode ierr;
103118c3aa7eSMark 
103218c3aa7eSMark   PetscFunctionBegin;
103318c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
103418c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
103518c3aa7eSMark   PetscFunctionReturn(0);
103618c3aa7eSMark }
103718c3aa7eSMark 
103818c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n)
103918c3aa7eSMark {
104018c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
104118c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
104218c3aa7eSMark 
104318c3aa7eSMark   PetscFunctionBegin;
104418c3aa7eSMark   pc_gamg->esteig_max_it = n;
104518c3aa7eSMark   PetscFunctionReturn(0);
104618c3aa7eSMark }
104718c3aa7eSMark 
104818c3aa7eSMark /*@
104918c3aa7eSMark    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother
105018c3aa7eSMark 
105118c3aa7eSMark    Collective on PC
105218c3aa7eSMark 
105318c3aa7eSMark    Input Parameters:
105418c3aa7eSMark +  pc - the preconditioner context
105518c3aa7eSMark -  n - number of its
105618c3aa7eSMark 
105718c3aa7eSMark    Options Database Key:
105818c3aa7eSMark .  -pc_gamg_use_sa_esteig <true,false>
105918c3aa7eSMark 
106018c3aa7eSMark    Notes:
106118c3aa7eSMark 
106218c3aa7eSMark    Level: intermediate
106318c3aa7eSMark 
106418c3aa7eSMark .seealso: ()
106518c3aa7eSMark @*/
106618c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
106718c3aa7eSMark {
106818c3aa7eSMark   PetscErrorCode ierr;
106918c3aa7eSMark 
107018c3aa7eSMark   PetscFunctionBegin;
107118c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
107218c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
107318c3aa7eSMark   PetscFunctionReturn(0);
107418c3aa7eSMark }
107518c3aa7eSMark 
10760ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
107718c3aa7eSMark {
107818c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
107918c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
108018c3aa7eSMark 
108118c3aa7eSMark   PetscFunctionBegin;
108218c3aa7eSMark   pc_gamg->use_sa_esteig = n ? 1 : 0;
108318c3aa7eSMark   PetscFunctionReturn(0);
108418c3aa7eSMark }
108518c3aa7eSMark 
108618c3aa7eSMark /*@C
108718c3aa7eSMark    PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby)
108818c3aa7eSMark 
108918c3aa7eSMark    Collective on PC
109018c3aa7eSMark 
109118c3aa7eSMark    Input Parameters:
109218c3aa7eSMark +  pc - the preconditioner context
109318c3aa7eSMark -  t - ksp type
109418c3aa7eSMark 
109518c3aa7eSMark    Options Database Key:
109618c3aa7eSMark .  -pc_gamg_esteig_ksp_type <type>
109718c3aa7eSMark 
109818c3aa7eSMark    Notes:
109918c3aa7eSMark 
110018c3aa7eSMark    Level: intermediate
110118c3aa7eSMark 
110218c3aa7eSMark .seealso: ()
110318c3aa7eSMark @*/
110418c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[])
110518c3aa7eSMark {
110618c3aa7eSMark   PetscErrorCode ierr;
110718c3aa7eSMark 
110818c3aa7eSMark   PetscFunctionBegin;
110918c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111018c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr);
111118c3aa7eSMark   PetscFunctionReturn(0);
111218c3aa7eSMark }
111318c3aa7eSMark 
111418c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[])
111518c3aa7eSMark {
111618c3aa7eSMark   PetscErrorCode ierr;
111718c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
111818c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
111918c3aa7eSMark 
112018c3aa7eSMark   PetscFunctionBegin;
112118c3aa7eSMark   ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr);
112218c3aa7eSMark   PetscFunctionReturn(0);
112318c3aa7eSMark }
112418c3aa7eSMark 
112518c3aa7eSMark /*@
112618c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
112718c3aa7eSMark 
112818c3aa7eSMark    Collective on PC
112918c3aa7eSMark 
113018c3aa7eSMark    Input Parameters:
113118c3aa7eSMark +  pc - the preconditioner context
113218c3aa7eSMark -  emax - max eigenvalue
113318c3aa7eSMark .  emin - min eigenvalue
113418c3aa7eSMark 
113518c3aa7eSMark    Options Database Key:
113618c3aa7eSMark .  -pc_gamg_eigenvalues
113718c3aa7eSMark 
113818c3aa7eSMark    Level: intermediate
113918c3aa7eSMark 
114018c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType()
114118c3aa7eSMark @*/
114218c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
114318c3aa7eSMark {
114418c3aa7eSMark   PetscErrorCode ierr;
114518c3aa7eSMark 
114618c3aa7eSMark   PetscFunctionBegin;
114718c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
114818c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr);
114918c3aa7eSMark   PetscFunctionReturn(0);
115018c3aa7eSMark }
115118c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
115218c3aa7eSMark {
115318c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
115418c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
115518c3aa7eSMark 
115618c3aa7eSMark   PetscFunctionBegin;
115718c3aa7eSMark   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);
115818c3aa7eSMark   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);
115918c3aa7eSMark   pc_gamg->emax = emax;
116018c3aa7eSMark   pc_gamg->emin = emin;
116118c3aa7eSMark 
116218c3aa7eSMark   PetscFunctionReturn(0);
116318c3aa7eSMark }
116418c3aa7eSMark 
116518c3aa7eSMark /*@
1166cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams    Collective on PC
1169dfd5c07aSMark F. Adams 
1170dfd5c07aSMark F. Adams    Input Parameters:
11711cc46a46SBarry Smith +  pc - the preconditioner context
11721cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1173dfd5c07aSMark F. Adams 
1174dfd5c07aSMark F. Adams    Options Database Key:
11751cc46a46SBarry Smith .  -pc_gamg_reuse_interpolation <true,false>
1176dfd5c07aSMark F. Adams 
1177dfd5c07aSMark F. Adams    Level: intermediate
1178dfd5c07aSMark F. Adams 
117995452b02SPatrick Sanan    Notes:
118095452b02SPatrick Sanan     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1181cab9ed1eSBarry Smith           rebuilding the preconditioner quicker.
1182cab9ed1eSBarry Smith 
1183dfd5c07aSMark F. Adams .seealso: ()
1184dfd5c07aSMark F. Adams @*/
11851cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1186dfd5c07aSMark F. Adams {
1187dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1188dfd5c07aSMark F. Adams 
1189dfd5c07aSMark F. Adams   PetscFunctionBegin;
1190dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11911cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1192dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1193dfd5c07aSMark F. Adams }
1194dfd5c07aSMark F. Adams 
11951cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1196dfd5c07aSMark F. Adams {
1197dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1198dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1199dfd5c07aSMark F. Adams 
1200dfd5c07aSMark F. Adams   PetscFunctionBegin;
1201dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1202dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1203dfd5c07aSMark F. Adams }
1204dfd5c07aSMark F. Adams 
1205ffc955d6SMark F. Adams /*@
1206cab9ed1eSBarry 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.
1207ffc955d6SMark F. Adams 
1208ffc955d6SMark F. Adams    Collective on PC
1209ffc955d6SMark F. Adams 
1210ffc955d6SMark F. Adams    Input Parameters:
1211cab9ed1eSBarry Smith +  pc - the preconditioner context
1212cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1213ffc955d6SMark F. Adams 
1214ffc955d6SMark F. Adams    Options Database Key:
1215cab9ed1eSBarry Smith .  -pc_gamg_asm_use_agg
1216ffc955d6SMark F. Adams 
1217ffc955d6SMark F. Adams    Level: intermediate
1218ffc955d6SMark F. Adams 
1219ffc955d6SMark F. Adams .seealso: ()
1220ffc955d6SMark F. Adams @*/
1221cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1222ffc955d6SMark F. Adams {
1223ffc955d6SMark F. Adams   PetscErrorCode ierr;
1224ffc955d6SMark F. Adams 
1225ffc955d6SMark F. Adams   PetscFunctionBegin;
1226ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1227cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1228ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1229ffc955d6SMark F. Adams }
1230ffc955d6SMark F. Adams 
1231cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1232ffc955d6SMark F. Adams {
1233ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1234ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1235ffc955d6SMark F. Adams 
1236ffc955d6SMark F. Adams   PetscFunctionBegin;
1237cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1238ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1239ffc955d6SMark F. Adams }
1240ffc955d6SMark F. Adams 
1241171cca9aSMark Adams /*@
1242cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1243171cca9aSMark Adams 
1244171cca9aSMark Adams    Collective on PC
1245171cca9aSMark Adams 
1246171cca9aSMark Adams    Input Parameters:
1247171cca9aSMark Adams +  pc - the preconditioner context
1248cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1249171cca9aSMark Adams 
1250171cca9aSMark Adams    Options Database Key:
1251cf8ae1d3SMark Adams .  -pc_gamg_use_parallel_coarse_grid_solver
1252171cca9aSMark Adams 
1253171cca9aSMark Adams    Level: intermediate
1254171cca9aSMark Adams 
125539d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1256171cca9aSMark Adams @*/
1257171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1258171cca9aSMark Adams {
1259171cca9aSMark Adams   PetscErrorCode ierr;
1260171cca9aSMark Adams 
1261171cca9aSMark Adams   PetscFunctionBegin;
1262171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1263171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1264171cca9aSMark Adams   PetscFunctionReturn(0);
1265171cca9aSMark Adams }
1266171cca9aSMark Adams 
1267171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1268171cca9aSMark Adams {
1269171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1270171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1271171cca9aSMark Adams 
1272171cca9aSMark Adams   PetscFunctionBegin;
1273171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1274ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1275ffc955d6SMark F. Adams }
1276ffc955d6SMark F. Adams 
12774ef23d27SMark F. Adams /*@
1278ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1279ce7c7f2fSMark Adams 
1280ce7c7f2fSMark Adams    Collective on PC
1281ce7c7f2fSMark Adams 
1282ce7c7f2fSMark Adams    Input Parameters:
1283ce7c7f2fSMark Adams +  pc - the preconditioner context
1284ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1285ce7c7f2fSMark Adams 
1286ce7c7f2fSMark Adams    Options Database Key:
1287ce7c7f2fSMark Adams .  -pc_gamg_cpu_pin_coarse_grids
1288ce7c7f2fSMark Adams 
1289ce7c7f2fSMark Adams    Level: intermediate
1290ce7c7f2fSMark Adams 
129139d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1292ce7c7f2fSMark Adams @*/
1293ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1294ce7c7f2fSMark Adams {
1295ce7c7f2fSMark Adams   PetscErrorCode ierr;
1296ce7c7f2fSMark Adams 
1297ce7c7f2fSMark Adams   PetscFunctionBegin;
1298ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1299ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1300ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1301ce7c7f2fSMark Adams }
1302ce7c7f2fSMark Adams 
1303ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1304ce7c7f2fSMark Adams {
1305ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1306ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1307ce7c7f2fSMark Adams 
1308ce7c7f2fSMark Adams   PetscFunctionBegin;
1309ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1310ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1311ce7c7f2fSMark Adams }
1312ce7c7f2fSMark Adams 
1313ce7c7f2fSMark Adams /*@
1314ce7c7f2fSMark Adams    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1315ce7c7f2fSMark Adams 
1316ce7c7f2fSMark Adams    Collective on PC
1317ce7c7f2fSMark Adams 
1318ce7c7f2fSMark Adams    Input Parameters:
1319ce7c7f2fSMark Adams +  pc - the preconditioner context
1320ce7c7f2fSMark Adams -  flg - Layout type
1321ce7c7f2fSMark Adams 
1322ce7c7f2fSMark Adams    Options Database Key:
1323ce7c7f2fSMark Adams .  -pc_gamg_coarse_grid_layout_type
1324ce7c7f2fSMark Adams 
1325ce7c7f2fSMark Adams    Level: intermediate
1326ce7c7f2fSMark Adams 
132739d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1328ce7c7f2fSMark Adams @*/
1329ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1330ce7c7f2fSMark Adams {
1331ce7c7f2fSMark Adams   PetscErrorCode ierr;
1332ce7c7f2fSMark Adams 
1333ce7c7f2fSMark Adams   PetscFunctionBegin;
1334ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1335ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1336ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1337ce7c7f2fSMark Adams }
1338ce7c7f2fSMark Adams 
1339ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1340ce7c7f2fSMark Adams {
1341ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1342ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1343ce7c7f2fSMark Adams 
1344ce7c7f2fSMark Adams   PetscFunctionBegin;
1345ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1346ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1347ce7c7f2fSMark Adams }
1348ce7c7f2fSMark Adams 
1349ce7c7f2fSMark Adams /*@
13501cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
13514ef23d27SMark F. Adams 
13524ef23d27SMark F. Adams    Not collective on PC
13534ef23d27SMark F. Adams 
13544ef23d27SMark F. Adams    Input Parameters:
13551cc46a46SBarry Smith +  pc - the preconditioner
13561cc46a46SBarry Smith -  n - the maximum number of levels to use
13574ef23d27SMark F. Adams 
13584ef23d27SMark F. Adams    Options Database Key:
13594ef23d27SMark F. Adams .  -pc_mg_levels
13604ef23d27SMark F. Adams 
13614ef23d27SMark F. Adams    Level: intermediate
13624ef23d27SMark F. Adams 
13634ef23d27SMark F. Adams .seealso: ()
13644ef23d27SMark F. Adams @*/
13654ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
13664ef23d27SMark F. Adams {
13674ef23d27SMark F. Adams   PetscErrorCode ierr;
13684ef23d27SMark F. Adams 
13694ef23d27SMark F. Adams   PetscFunctionBegin;
13704ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13714ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
13724ef23d27SMark F. Adams   PetscFunctionReturn(0);
13734ef23d27SMark F. Adams }
13744ef23d27SMark F. Adams 
13751e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
13764ef23d27SMark F. Adams {
13774ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
13784ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
13794ef23d27SMark F. Adams 
13804ef23d27SMark F. Adams   PetscFunctionBegin;
13819d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
13824ef23d27SMark F. Adams   PetscFunctionReturn(0);
13834ef23d27SMark F. Adams }
13844ef23d27SMark F. Adams 
13853542efc5SMark F. Adams /*@
13863542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
13873542efc5SMark F. Adams 
13883542efc5SMark F. Adams    Not collective on PC
13893542efc5SMark F. Adams 
13903542efc5SMark F. Adams    Input Parameters:
13911cc46a46SBarry Smith +  pc - the preconditioner context
1392c9567895SMark .  v - 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
1393055c8bd0SJed Brown -  n - number of threshold values provided in array
13943542efc5SMark F. Adams 
13953542efc5SMark F. Adams    Options Database Key:
13961cc46a46SBarry Smith .  -pc_gamg_threshold <threshold>
13973542efc5SMark F. Adams 
139895452b02SPatrick Sanan    Notes:
1399af3c827dSMark 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.
1400af3c827dSMark 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.
1401cab9ed1eSBarry Smith 
1402055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1403055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1404055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1405055c8bd0SJed Brown 
14063542efc5SMark F. Adams    Level: intermediate
14073542efc5SMark F. Adams 
1408af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
14093542efc5SMark F. Adams @*/
1410c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
14113542efc5SMark F. Adams {
14123542efc5SMark F. Adams   PetscErrorCode ierr;
14133542efc5SMark F. Adams 
14143542efc5SMark F. Adams   PetscFunctionBegin;
14153542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1416055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
1417c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
14183542efc5SMark F. Adams   PetscFunctionReturn(0);
14193542efc5SMark F. Adams }
14203542efc5SMark F. Adams 
1421c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
14223542efc5SMark F. Adams {
1423c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1424c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1425c1eae691SMark Adams   PetscInt i;
1426c1eae691SMark Adams   PetscFunctionBegin;
1427055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1428055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1429c1eae691SMark Adams   PetscFunctionReturn(0);
1430c1eae691SMark Adams }
1431c1eae691SMark Adams 
1432c1eae691SMark Adams /*@
1433c9567895SMark    PCGAMGSetRankReductionFactors - Set manual schedual for process reduction on coarse grids
1434c9567895SMark 
1435c9567895SMark    Collective on PC
1436c9567895SMark 
1437c9567895SMark    Input Parameters:
1438c9567895SMark +  pc - the preconditioner context
1439c9567895SMark .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1440c9567895SMark -  n - number of values provided in array
1441c9567895SMark 
1442c9567895SMark    Options Database Key:
1443c9567895SMark .  -pc_gamg_rank_reduction_factors <factors>
1444c9567895SMark 
1445c9567895SMark    Level: intermediate
1446c9567895SMark 
1447c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1448c9567895SMark @*/
1449c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1450c9567895SMark {
1451c9567895SMark   PetscErrorCode ierr;
1452c9567895SMark 
1453c9567895SMark   PetscFunctionBegin;
1454c9567895SMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1455c9567895SMark   if (n) PetscValidIntPointer(v,2);
1456c9567895SMark   ierr = PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));CHKERRQ(ierr);
1457c9567895SMark   PetscFunctionReturn(0);
1458c9567895SMark }
1459c9567895SMark 
1460c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1461c9567895SMark {
1462c9567895SMark   PC_MG   *mg      = (PC_MG*)pc->data;
1463c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1464c9567895SMark   PetscInt i;
1465c9567895SMark   PetscFunctionBegin;
1466c9567895SMark   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1467c9567895SMark   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1468c9567895SMark   PetscFunctionReturn(0);
1469c9567895SMark }
1470c9567895SMark 
1471c9567895SMark /*@
1472c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1473c1eae691SMark Adams 
1474c1eae691SMark Adams    Not collective on PC
1475c1eae691SMark Adams 
1476c1eae691SMark Adams    Input Parameters:
1477c1eae691SMark Adams +  pc - the preconditioner context
1478c1eae691SMark Adams -  scale - the threshold value reduction, ussually < 1.0
1479c1eae691SMark Adams 
1480c1eae691SMark Adams    Options Database Key:
1481c1eae691SMark Adams .  -pc_gamg_threshold_scale <v>
1482c1eae691SMark Adams 
1483055c8bd0SJed Brown    Notes:
1484055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1485055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1486055c8bd0SJed Brown 
1487c1eae691SMark Adams    Level: advanced
1488c1eae691SMark Adams 
1489055c8bd0SJed Brown .seealso: PCGAMGSetThreshold()
1490c1eae691SMark Adams @*/
1491c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1492c1eae691SMark Adams {
1493c1eae691SMark Adams   PetscErrorCode ierr;
14943542efc5SMark F. Adams 
14953542efc5SMark F. Adams   PetscFunctionBegin;
1496c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1497c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1498c1eae691SMark Adams   PetscFunctionReturn(0);
1499c1eae691SMark Adams }
1500c1eae691SMark Adams 
1501c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1502c1eae691SMark Adams {
1503c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1504c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1505c1eae691SMark Adams   PetscFunctionBegin;
1506c1eae691SMark Adams   pc_gamg->threshold_scale = v;
15073542efc5SMark F. Adams   PetscFunctionReturn(0);
15083542efc5SMark F. Adams }
15093542efc5SMark F. Adams 
1510e20c40e8SBarry Smith /*@C
1511c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1512676e1743SMark F. Adams 
1513676e1743SMark F. Adams    Collective on PC
1514676e1743SMark F. Adams 
1515676e1743SMark F. Adams    Input Parameters:
1516c60c7ad4SBarry Smith +  pc - the preconditioner context
1517c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1518676e1743SMark F. Adams 
1519676e1743SMark F. Adams    Options Database Key:
1520cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1521676e1743SMark F. Adams 
1522676e1743SMark F. Adams    Level: intermediate
1523676e1743SMark F. Adams 
1524cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1525676e1743SMark F. Adams @*/
152619fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1527676e1743SMark F. Adams {
1528676e1743SMark F. Adams   PetscErrorCode ierr;
1529676e1743SMark F. Adams 
1530676e1743SMark F. Adams   PetscFunctionBegin;
1531676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1532806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1533676e1743SMark F. Adams   PetscFunctionReturn(0);
1534676e1743SMark F. Adams }
1535676e1743SMark F. Adams 
1536e20c40e8SBarry Smith /*@C
1537c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1538c60c7ad4SBarry Smith 
1539c60c7ad4SBarry Smith    Collective on PC
1540c60c7ad4SBarry Smith 
1541c60c7ad4SBarry Smith    Input Parameter:
1542c60c7ad4SBarry Smith .  pc - the preconditioner context
1543c60c7ad4SBarry Smith 
1544c60c7ad4SBarry Smith    Output Parameter:
1545c60c7ad4SBarry Smith .  type - the type of algorithm used
1546c60c7ad4SBarry Smith 
1547c60c7ad4SBarry Smith    Level: intermediate
1548c60c7ad4SBarry Smith 
15491c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1550c60c7ad4SBarry Smith @*/
1551c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1552c60c7ad4SBarry Smith {
1553c60c7ad4SBarry Smith   PetscErrorCode ierr;
1554c60c7ad4SBarry Smith 
1555c60c7ad4SBarry Smith   PetscFunctionBegin;
1556c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1557c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1558c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1559c60c7ad4SBarry Smith }
1560c60c7ad4SBarry Smith 
1561c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1562c60c7ad4SBarry Smith {
1563c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1564c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1565c60c7ad4SBarry Smith 
1566c60c7ad4SBarry Smith   PetscFunctionBegin;
1567c60c7ad4SBarry Smith   *type = pc_gamg->type;
1568c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1569c60c7ad4SBarry Smith }
1570c60c7ad4SBarry Smith 
15711e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1572676e1743SMark F. Adams {
15739d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
15741ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
15751ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1576676e1743SMark F. Adams 
1577676e1743SMark F. Adams   PetscFunctionBegin;
1578c60c7ad4SBarry Smith   pc_gamg->type = type;
15791c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
15809d5b6da9SMark F. Adams   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
15811ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
15821ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
15831ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1584e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
15853ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
15863ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
15873ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
15883ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
15893ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
15903ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
15913ae0bb68SMark Adams     pc_gamg->data_sz = 0;
15921ab5ffc9SJed Brown   }
15931ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
15941ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
15959d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1596676e1743SMark F. Adams   PetscFunctionReturn(0);
1597676e1743SMark F. Adams }
1598676e1743SMark F. Adams 
15999d766c59SMark Adams /* -------------------------------------------------------------------------- */
16009d766c59SMark Adams /*
16019d766c59SMark Adams    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
16029d766c59SMark Adams 
16039d766c59SMark Adams    Input Parameter:
16049d766c59SMark Adams .  pc - the preconditioner context
16059d766c59SMark Adams 
16069d766c59SMark Adams    Output Parameter:
16079d766c59SMark Adams .  gc - grid complexity = sum_i(nnz_i) / nnz_0
16089d766c59SMark Adams 
16099d766c59SMark Adams    Level: advanced
16109d766c59SMark Adams */
16119d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
16129d766c59SMark Adams {
16139d766c59SMark Adams   PetscErrorCode ierr;
16149d766c59SMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
16159d766c59SMark Adams   PC_MG_Levels   **mglevels = mg->levels;
16163966268fSBarry Smith   PetscInt       lev;
16173966268fSBarry Smith   PetscLogDouble nnz0 = 0, sgc = 0;
16189d766c59SMark Adams   MatInfo        info;
16193966268fSBarry Smith 
16209d766c59SMark Adams   PetscFunctionBegin;
1621dbf6bb8dSprj-   if (!pc->setupcalled) {
1622dbf6bb8dSprj-     *gc = 0;
1623dbf6bb8dSprj-     PetscFunctionReturn(0);
1624dbf6bb8dSprj-   }
16259d766c59SMark Adams   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
16263966268fSBarry Smith   for (lev=0; lev<mg->nlevels; lev++) {
162762a6e064SMark Adams     Mat dB;
162862a6e064SMark Adams     ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr);
162962a6e064SMark Adams     ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
16303966268fSBarry Smith     sgc += info.nz_used;
16319d766c59SMark Adams     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
16329d766c59SMark Adams   }
16333966268fSBarry Smith   if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0);
16343966268fSBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
16359d766c59SMark Adams   PetscFunctionReturn(0);
16369d766c59SMark Adams }
16379d766c59SMark Adams 
16385adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
16395adeb434SBarry Smith {
1640c1eae691SMark Adams   PetscErrorCode ierr,i;
16415adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
16425adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
164323b2d91dSMark Adams   PetscReal       gc=0;
16445adeb434SBarry Smith   PetscFunctionBegin;
16455adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1646459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1647b3e187dcSStefano Zampini   for (i=0;i<mg->nlevels; i++) {
1648c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1649c1eae691SMark Adams   }
1650459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1651459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1652cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1653cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1654cab9ed1eSBarry Smith   }
1655171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1656171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1657171cca9aSMark Adams   }
1658ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1659ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1660ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1661ce7c7f2fSMark Adams   }
1662ce7c7f2fSMark Adams #endif
1663ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1664ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1665ce7c7f2fSMark Adams   /* } else { */
1666ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1667ce7c7f2fSMark Adams   /* } */
16685adeb434SBarry Smith   if (pc_gamg->ops->view) {
16695adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
16705adeb434SBarry Smith   }
16719d766c59SMark Adams   ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr);
16729d766c59SMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);CHKERRQ(ierr);
16735adeb434SBarry Smith   PetscFunctionReturn(0);
16745adeb434SBarry Smith }
16755adeb434SBarry Smith 
16764416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
16775b89ad90SMark F. Adams {
1678676e1743SMark F. Adams   PetscErrorCode ierr;
1679676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1680676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
168118c3aa7eSMark   PetscBool      flag,f2;
16823b4367a7SBarry Smith   MPI_Comm       comm;
168318c3aa7eSMark   char           prefix[256],tname[32];
1684c1eae691SMark Adams   PetscInt       i,n;
168514a9496bSBarry Smith   const char     *pcpre;
16860a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
16875b89ad90SMark F. Adams   PetscFunctionBegin;
16883b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1689e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
16901a1c1e04SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1691bd94a7aaSJed Brown   if (flag) {
1692bd94a7aaSJed Brown     ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
16931ab5ffc9SJed Brown   }
169418c3aa7eSMark   ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr);
169518c3aa7eSMark   if (flag) {
169618c3aa7eSMark     ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr);
169718c3aa7eSMark   }
1698cab9ed1eSBarry Smith   ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
169918c3aa7eSMark   f2 = PETSC_TRUE;
170018c3aa7eSMark   ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr);
170118c3aa7eSMark   if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0;
17021cc46a46SBarry Smith   ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1703a303c832SJed 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);
1704cf8ae1d3SMark 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);
1705ce7c7f2fSMark 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);
1706a0095786SMark   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);
170794ae4db5SBarry 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);
170818c3aa7eSMark   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);
170994ae4db5SBarry 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);
1710a303c832SJed 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);
171118c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
1712c1eae691SMark Adams   ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
171318c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1714efd3c5ceSMark Adams     if (!flag) n = 1;
1715c1eae691SMark Adams     i = n;
171618c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1717c1eae691SMark Adams   }
1718c9567895SMark   n = PETSC_MG_MAXLEVELS;
1719c9567895SMark   ierr = PetscOptionsIntArray("-pc_gamg_rank_reduction_factors","Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)","PCGAMGSetRankReductionFactors",pc_gamg->level_reduction_factors,&n,&flag);CHKERRQ(ierr);
1720c9567895SMark   if (!flag) i = 0;
1721c9567895SMark   else i = n;
1722c9567895SMark   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
172394ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
172418c3aa7eSMark   {
172518c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
172618c3aa7eSMark     n = 2;
172718c3aa7eSMark     ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr);
172818c3aa7eSMark     if (flag) {
172918c3aa7eSMark       if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
173018c3aa7eSMark       ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr);
173118c3aa7eSMark     }
173218c3aa7eSMark   }
1733b7cbab4eSMark Adams   /* set options for subtype */
1734e55864a3SBarry Smith   if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
173518c3aa7eSMark 
173614a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
173714a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1738676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
17395b89ad90SMark F. Adams   PetscFunctionReturn(0);
17405b89ad90SMark F. Adams }
17415b89ad90SMark F. Adams 
17425b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
17435b89ad90SMark F. Adams /*MC
17441cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
17455b89ad90SMark F. Adams 
1746280d9858SJed Brown    Options Database Keys:
1747cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1748cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1749cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1750cab9ed1eSBarry 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
1751cab9ed1eSBarry 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>
1752cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1753cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
17546008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1755c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1756cab9ed1eSBarry Smith 
1757cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1758cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1759cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1760cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1761cab9ed1eSBarry Smith 
1762db9745e2SBarry Smith    Multigrid options:
1763db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1764db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1765db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1766cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
17675b89ad90SMark F. Adams 
17681cc46a46SBarry Smith 
176995452b02SPatrick Sanan   Notes:
177095452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1771db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1772db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1773db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
17741cc46a46SBarry Smith 
17755b89ad90SMark F. Adams   Level: intermediate
1776280d9858SJed Brown 
17771cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
177818c3aa7eSMark            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType()
17795b89ad90SMark F. Adams M*/
1780b2573a8aSBarry Smith 
17818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
17825b89ad90SMark F. Adams {
1783c1eae691SMark Adams   PetscErrorCode ierr,i;
17845b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
17855b89ad90SMark F. Adams   PC_MG          *mg;
17865b89ad90SMark F. Adams 
17875b89ad90SMark F. Adams   PetscFunctionBegin;
17881c1aac46SBarry Smith    /* register AMG type */
17891c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
17901c1aac46SBarry Smith 
17915b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
17921c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
17935b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
17945b89ad90SMark F. Adams 
17955b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1796b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
179769aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
17985b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
17995b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
18005b89ad90SMark F. Adams 
1801b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
18021ab5ffc9SJed Brown 
18039d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
18049d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
18050a545947SLisandro Dalcin   pc_gamg->data    = NULL;
18065b89ad90SMark F. Adams 
18079d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
18085b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
18095b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
18105b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
18115b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
18125adeb434SBarry Smith   mg->view                = PCView_GAMG;
18135b89ad90SMark F. Adams 
181497d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
181597d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1816bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1817bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1818cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
181918c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr);
182018c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr);
182118c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr);
182218c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr);
18231cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1824cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1825171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1826ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1827ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1828bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1829c9567895SMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);CHKERRQ(ierr);
1830c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1831bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1832c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1833bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
18349d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1835d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
18360c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1837171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1838a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1839a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1840038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
184125a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
184218c3aa7eSMark   for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1843c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
184418c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
18459ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
1846d24ecf33SMark   ierr = PetscStrcpy(pc_gamg->esteig_type,NULL);CHKERRQ(ierr);
184718c3aa7eSMark   pc_gamg->esteig_max_it    = 10;
184818c3aa7eSMark   pc_gamg->use_sa_esteig    = -1;
184918c3aa7eSMark   pc_gamg->emin             = 0;
185018c3aa7eSMark   pc_gamg->emax             = 0;
185118c3aa7eSMark 
1852c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
18539d5b6da9SMark F. Adams 
1854bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1855bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
18565b89ad90SMark F. Adams   PetscFunctionReturn(0);
18575b89ad90SMark F. Adams }
18583e3471ccSMark Adams 
18593e3471ccSMark Adams /*@C
18603e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
18618a690491SBarry Smith     from PCInitializePackage().
18623e3471ccSMark Adams 
18633e3471ccSMark Adams  Level: developer
18643e3471ccSMark Adams 
18653e3471ccSMark Adams  .seealso: PetscInitialize()
18663e3471ccSMark Adams @*/
18673e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
18683e3471ccSMark Adams {
18693e3471ccSMark Adams   PetscErrorCode ierr;
18703e3471ccSMark Adams 
18713e3471ccSMark Adams   PetscFunctionBegin;
18723e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
18733e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
18743e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
18753e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
18768e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
18773e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1878c1c463dbSMark Adams 
1879c1c463dbSMark Adams   /* general events */
1880fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1881fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1882fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1883fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1884c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1885c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1886fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1887c1c463dbSMark Adams 
18885b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG
18895b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
18905b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
18915b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
18925b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
18935b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
18945b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
18955b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
18965b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1897bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
18985b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
18995b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
19005b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
19015b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
19025b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
19035b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
19045b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
19055b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
19065b89ad90SMark F. Adams 
19075b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
19085b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
19095b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
19105b89ad90SMark F. Adams   /* create timer stages */
19115b89ad90SMark F. Adams #if defined GAMG_STAGES
19125b89ad90SMark F. Adams   {
19135b89ad90SMark F. Adams     char     str[32];
19145b89ad90SMark F. Adams     PetscInt lidx;
19155b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
19165b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
19175b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
1918c9567895SMark       sprintf(str,"MG Level %d",(int)lidx);
19195b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
19205b89ad90SMark F. Adams     }
19215b89ad90SMark F. Adams   }
19225b89ad90SMark F. Adams #endif
19235b89ad90SMark F. Adams #endif
19243e3471ccSMark Adams   PetscFunctionReturn(0);
19253e3471ccSMark Adams }
19263e3471ccSMark Adams 
19273e3471ccSMark Adams /*@C
19281c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
19291c1aac46SBarry Smith     called from PetscFinalize() automatically.
19303e3471ccSMark Adams 
19313e3471ccSMark Adams  Level: developer
19323e3471ccSMark Adams 
19333e3471ccSMark Adams  .seealso: PetscFinalize()
19343e3471ccSMark Adams @*/
19353e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
19363e3471ccSMark Adams {
19373e3471ccSMark Adams   PetscErrorCode ierr;
19383e3471ccSMark Adams 
19393e3471ccSMark Adams   PetscFunctionBegin;
19403e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
19413e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
19423e3471ccSMark Adams   PetscFunctionReturn(0);
19433e3471ccSMark Adams }
1944a36cf38bSToby Isaac 
1945a36cf38bSToby Isaac /*@C
1946a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1947a36cf38bSToby Isaac 
1948a36cf38bSToby Isaac  Input Parameters:
1949a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1950a36cf38bSToby Isaac  - create - function for creating the gamg context.
1951a36cf38bSToby Isaac 
1952a36cf38bSToby Isaac   Level: advanced
1953a36cf38bSToby Isaac 
19541c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1955a36cf38bSToby Isaac @*/
1956a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1957a36cf38bSToby Isaac {
1958a36cf38bSToby Isaac   PetscErrorCode ierr;
1959a36cf38bSToby Isaac 
1960a36cf38bSToby Isaac   PetscFunctionBegin;
1961a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1962a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1963a36cf38bSToby Isaac   PetscFunctionReturn(0);
1964a36cf38bSToby Isaac }
1965a36cf38bSToby Isaac 
1966