xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 147403d98689287189d69e47992b7c8152b2c9da)
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*/
618c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>    /*I "petscksp.h" I*/
7f96513f1SMatthew G Knepley 
8c9567895SMark #if defined(PETSC_HAVE_CUDA)
9c9567895SMark   #include <cuda_runtime.h>
10c9567895SMark #endif
11c9567895SMark 
12c9567895SMark #if defined(PETSC_HAVE_HIP)
13c9567895SMark   #include <hip/hip_runtime.h>
14c9567895SMark #endif
15c9567895SMark 
160cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET];
174555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
18fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG;
19fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO;
200cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG;
210cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO;
220cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG;
230cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO;
24fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG;
250cbbd2e1SMark F. Adams 
26b8fd24d8SMark F. Adams /* #define GAMG_STAGES */
274555aa8cSStefano Zampini #if defined(GAMG_STAGES)
2818c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
29b4fbaa2aSMark F. Adams #endif
30f96513f1SMatthew G Knepley 
310a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
323e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized;
339d5b6da9SMark F. Adams 
34d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */
35d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc)
36d3d6bff4SMark F. Adams {
3718c3aa7eSMark   PetscErrorCode ierr, level;
38d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
39d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
40d3d6bff4SMark F. Adams 
41d3d6bff4SMark F. Adams   PetscFunctionBegin;
4222a233eaSStefano Zampini   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
431c1aac46SBarry Smith   pc_gamg->data_sz = 0;
44878e152fSMark F. Adams   ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr);
4518c3aa7eSMark   for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
4618c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
4718c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
4818c3aa7eSMark   }
4918c3aa7eSMark   pc_gamg->emin = 0;
5018c3aa7eSMark   pc_gamg->emax = 0;
51a2f3521dSMark F. Adams   PetscFunctionReturn(0);
52a2f3521dSMark F. Adams }
53a2f3521dSMark F. Adams 
545b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
555b89ad90SMark F. Adams /*
56c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
57a147abb0SMark F. Adams      of active processors.
585b89ad90SMark F. Adams 
595b89ad90SMark F. Adams    Input Parameter:
60a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
61a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
629d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
63c5bfad50SMark F. Adams    . cr_bs - coarse block size
643530afc2SMark F. Adams    In/Output Parameter:
65a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
66afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6711e60469SMark F. Adams    Output Parameter:
683530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
695b89ad90SMark F. Adams */
705cb416c2SMark F. Adams 
71171cca9aSMark 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)
725b89ad90SMark F. Adams {
73a2f3521dSMark F. Adams   PetscErrorCode  ierr;
749d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
75486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
76a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
773b4367a7SBarry Smith   MPI_Comm        comm;
78c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
793ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
805b89ad90SMark F. Adams 
815b89ad90SMark F. Adams   PetscFunctionBegin;
823b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr);
83ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
84ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
85c5bfad50SMark F. Adams   ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr);
864555aa8cSStefano Zampini   ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr);
879d5b6da9SMark F. Adams   ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr);
884555aa8cSStefano Zampini   ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);CHKERRQ(ierr);
89038e3b61SMark F. Adams 
90ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
91ce7c7f2fSMark Adams 
923ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
930298fd71SBarry Smith   ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr);
943ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
953ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9673911c69SBarry Smith   } else {
973ae0bb68SMark Adams     PetscInt  bs;
983ae0bb68SMark Adams     ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr);
993ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
1003ae0bb68SMark Adams   }
101c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
102c9567895SMark   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. */
103c9567895SMark #if defined(PETSC_HAVE_CUDA)
104c9567895SMark     PetscShmComm pshmcomm;
105c9567895SMark     PetscMPIInt  locrank;
106c9567895SMark     MPI_Comm     loccomm;
107c9567895SMark     PetscInt     s_nnodes,r_nnodes, new_new_size;
108c9567895SMark     cudaError_t  cerr;
109c9567895SMark     int          devCount;
110c9567895SMark     ierr = PetscShmCommGet(comm,&pshmcomm);CHKERRQ(ierr);
111c9567895SMark     ierr = PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);CHKERRQ(ierr);
11255b25c41SPierre Jolivet     ierr = MPI_Comm_rank(loccomm, &locrank);CHKERRMPI(ierr);
113c9567895SMark     s_nnodes = !locrank;
11455b25c41SPierre Jolivet     ierr = MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr);
1152c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size%r_nnodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes);
116c9567895SMark     devCount = 0;
117c9567895SMark     cerr = cudaGetDeviceCount(&devCount);
118c9567895SMark     cudaGetLastError(); /* Reset the last error */
119c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
120c9567895SMark       new_new_size = r_nnodes * devCount;
121c9567895SMark       new_size = new_new_size;
1221bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes);CHKERRQ(ierr);
123c9567895SMark     } else {
1241bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: With Cuda but no device. Use heuristics.");CHKERRQ(ierr);
125c9567895SMark       goto HEURISTIC;
126c9567895SMark     }
127c9567895SMark #else
128c9567895SMark     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
129c9567895SMark #endif
130c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
1312c71b3e2SJacob Faibussowitsch     PetscCheckFalse(nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level],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]);
132c9567895SMark     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
1331bde8fbeSMark Adams     ierr = PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%D)\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);CHKERRQ(ierr);
134c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
135c9567895SMark     new_size = 1;
1361bde8fbeSMark Adams     ierr = PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size);CHKERRQ(ierr);
137c9567895SMark   } else {
138472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
139c9567895SMark #if defined(PETSC_HAVE_CUDA)
140c9567895SMark     HEURISTIC:
141c9567895SMark #endif
1420298fd71SBarry Smith     ierr     = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr);
143a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1447f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
145c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
1461bde8fbeSMark Adams     ierr = PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size);CHKERRQ(ierr);
147a2f3521dSMark F. Adams   }
1482e3501ffSMark Adams   if (new_size==nactive) {
149ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
150ce7c7f2fSMark Adams     if (new_size < size) {
151ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
1521bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive);CHKERRQ(ierr);
153ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
154b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr);
155b470e4b4SRichard Tran Mills         ierr = MatBindToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr);
156ce7c7f2fSMark Adams       }
157ce7c7f2fSMark Adams     }
158ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1592e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
160192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
161885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
16271959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
1632c71b3e2SJacob Faibussowitsch     PetscCheckFalse(ncrs_eq % cr_bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
164ce7c7f2fSMark Adams     /* get new_size and rfactor */
165ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
166ce7c7f2fSMark Adams       /* find factor */
167ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
168ce7c7f2fSMark Adams       else {
169ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
170ce7c7f2fSMark Adams         jj = -1;
171ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
172ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
173ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
174ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
175ce7c7f2fSMark Adams             if (fact > best_fact) {
176ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
177ce7c7f2fSMark Adams             }
178ce7c7f2fSMark Adams           }
179ce7c7f2fSMark Adams         }
180ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
181ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
182ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
183ce7c7f2fSMark Adams         else expand_factor = rfactor;
184ce7c7f2fSMark Adams       }
185ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
1864cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1874cdfd227SMark         *a_Amat_crs = Cmat;
1881bde8fbeSMark Adams         ierr = PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",((PetscObject)pc)->prefix,new_size,ncrs_eq);CHKERRQ(ierr);
189ce7c7f2fSMark Adams         PetscFunctionReturn(0);
190ce7c7f2fSMark Adams       }
191ce7c7f2fSMark Adams     }
1924cdfd227SMark     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
193a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
194785e854fSJed Brown     ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr);
1952e3501ffSMark Adams     if (pc_gamg->repart) {
196a5b23f4aSJose E. Roman       /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
1975a9b9e01SMark F. Adams       Mat      adj;
1981bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",((PetscObject)pc)->prefix,*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");CHKERRQ(ierr);
199a2f3521dSMark F. Adams       /* get 'adj' */
200c5bfad50SMark F. Adams       if (cr_bs == 1) {
201038e3b61SMark F. Adams         ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
202806fa848SBarry Smith       } else {
203a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
204eb07cef2SMark F. Adams         Mat               tMat;
205a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
206b4fbaa2aSMark F. Adams         const PetscScalar *vals;
207b4fbaa2aSMark F. Adams         const PetscInt    *idx;
208a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
20939d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
210d9558ea9SBarry Smith         MatType           mtype;
211b4fbaa2aSMark F. Adams 
212e632b94dSBarry Smith         ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr);
213a2f3521dSMark F. Adams         ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr);
214a2f3521dSMark F. Adams         ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr);
215c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
2160a545947SLisandro Dalcin           ierr      = MatGetRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
217c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
218c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
2190a545947SLisandro Dalcin           ierr      = MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
2203ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
2213ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
22258471d46SMark F. Adams         }
2236876a03eSMark F. Adams 
224d9558ea9SBarry Smith         ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr);
2253b4367a7SBarry Smith         ierr = MatCreate(comm, &tMat);CHKERRQ(ierr);
2263ae0bb68SMark Adams         ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
227d9558ea9SBarry Smith         ierr = MatSetType(tMat,mtype);CHKERRQ(ierr);
228a2f3521dSMark F. Adams         ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr);
229a2f3521dSMark F. Adams         ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
230e632b94dSBarry Smith         ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
231eb07cef2SMark F. Adams 
232a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
233c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
23422063be5SMark F. Adams           ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
235eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
236c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
237eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
238eb07cef2SMark F. Adams             ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr);
239eb07cef2SMark F. Adams           }
24022063be5SMark F. Adams           ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr);
241eb07cef2SMark F. Adams         }
242eb07cef2SMark F. Adams         ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243eb07cef2SMark F. Adams         ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244eb07cef2SMark F. Adams 
245b4fbaa2aSMark F. Adams         if (llev++ == -1) {
246b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2478caf3d72SBarry Smith           ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr);
2483b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
249b4fbaa2aSMark F. Adams           ierr = MatView(tMat, viewer);CHKERRQ(ierr);
2503bf036e2SBarry Smith           ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
251b4fbaa2aSMark F. Adams         }
252eb07cef2SMark F. Adams         ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr);
253eb07cef2SMark F. Adams         ierr = MatDestroy(&tMat);CHKERRQ(ierr);
254a2f3521dSMark F. Adams       } /* create 'adj' */
255f150b916SMark F. Adams 
256a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2575a9b9e01SMark F. Adams         char            prefix[256];
2585a9b9e01SMark F. Adams         const char      *pcpre;
259b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
260b4fbaa2aSMark F. Adams         MatPartitioning mpart;
261a4b7d37bSMark F. Adams         IS              proc_is;
2622f03bc48SMark F. Adams 
2633b4367a7SBarry Smith         ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr);
2645ef31b24SMark F. Adams         ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr);
2659d5b6da9SMark F. Adams         ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
2668caf3d72SBarry Smith         ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
26759a0be82SJed Brown         ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
26811e60469SMark F. Adams         ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
269c5df96a5SBarry Smith         ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr);
270a4b7d37bSMark F. Adams         ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr);
27111e60469SMark F. Adams         ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
2725a9b9e01SMark F. Adams 
2735ef31b24SMark F. Adams         /* collect IS info */
274785e854fSJed Brown         ierr     = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr);
275a4b7d37bSMark F. Adams         ierr     = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr);
276a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
277c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
278ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
279eb07cef2SMark F. Adams           }
2805ef31b24SMark F. Adams         }
281a4b7d37bSMark F. Adams         ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr);
282a4b7d37bSMark F. Adams         ierr = ISDestroy(&proc_is);CHKERRQ(ierr);
2835ef31b24SMark F. Adams       }
2845ef31b24SMark F. Adams       ierr = MatDestroy(&adj);CHKERRQ(ierr);
2855a9b9e01SMark F. Adams 
2863b4367a7SBarry Smith       ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr);
2878263b398SMark F. Adams       ierr = PetscFree(newproc_idx);CHKERRQ(ierr);
28831cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
289ce7c7f2fSMark Adams       PetscInt targetPE;
2902c71b3e2SJacob Faibussowitsch       PetscCheckFalse(new_size==nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
2911bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: Number of equations (loc) %D with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq);CHKERRQ(ierr);
292ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2933b4367a7SBarry Smith       ierr     = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr);
294a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
295e33ef3b1SMark F. Adams 
29611e60469SMark F. Adams     /*
297a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
29811e60469SMark F. Adams     */
299a2f3521dSMark F. Adams     ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr);
3007700e67bSMark Adams     is_eq_num_prim = is_eq_num;
30111e60469SMark F. Adams     /*
302a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
30311e60469SMark F. Adams     */
304c5df96a5SBarry Smith     ierr        = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr);
305c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
306a2f3521dSMark F. Adams     ierr        = ISDestroy(&is_eq_newproc);CHKERRQ(ierr);
307ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
308a2f3521dSMark F. Adams 
309a2f3521dSMark F. Adams     ierr = PetscFree(counts);CHKERRQ(ierr);
310885364a3SMark 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 */
311885364a3SMark Adams     {
312885364a3SMark Adams       Vec            src_crd, dest_crd;
313885364a3SMark 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;
314885364a3SMark Adams       VecScatter     vecscat;
315885364a3SMark Adams       PetscScalar    *array;
316885364a3SMark Adams       IS isscat;
317a2f3521dSMark F. Adams       /* move data (for primal equations only) */
31822063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
3193b4367a7SBarry Smith       ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr);
3203ae0bb68SMark Adams       ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr);
321c0dedaeaSBarry Smith       ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */
32211e60469SMark F. Adams       /*
3239d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
324c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
32511e60469SMark F. Adams       */
326854ce69bSBarry Smith       ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr);
327a2f3521dSMark F. Adams       ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3283ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
329c5bfad50SMark F. Adams         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
330a2f3521dSMark F. Adams         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
33111e60469SMark F. Adams       }
332a2f3521dSMark F. Adams       ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr);
3333ae0bb68SMark Adams       ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr);
33492a756f0SMark F. Adams       ierr = PetscFree(tidx);CHKERRQ(ierr);
33511e60469SMark F. Adams       /*
33611e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
33711e60469SMark F. Adams       */
3383ae0bb68SMark Adams       ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr);
3399d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3403ae0bb68SMark Adams         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3413ae0bb68SMark Adams         for (ii=0; ii<ncrs; ii++) {
3429d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
343a2f3521dSMark F. Adams             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
344c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
345676e1743SMark F. Adams             ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr);
346d3d6bff4SMark F. Adams           }
347038e3b61SMark F. Adams         }
348eb07cef2SMark F. Adams       }
349eb07cef2SMark F. Adams       ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr);
350eb07cef2SMark F. Adams       ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr);
35111e60469SMark F. Adams       /*
35211e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
35311e60469SMark F. Adams         to the correct processor
35411e60469SMark F. Adams       */
3559448b7f1SJunchao Zhang       ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr);
35611e60469SMark F. Adams       ierr = ISDestroy(&isscat);CHKERRQ(ierr);
35711e60469SMark F. Adams       ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
35811e60469SMark F. Adams       ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
35911e60469SMark F. Adams       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
36011e60469SMark F. Adams       ierr = VecDestroy(&src_crd);CHKERRQ(ierr);
36111e60469SMark F. Adams       /*
36211e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
36311e60469SMark F. Adams       */
364c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
365578f55a3SPeter Brune       ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr);
3662fa5cd67SKarl Rupp 
3673ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3683ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3692fa5cd67SKarl Rupp 
37011e60469SMark F. Adams       ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr);
3719d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3723ae0bb68SMark Adams         for (ii=0; ii<ncrs_new; ii++) {
3739d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
374a2f3521dSMark F. Adams             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
375c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
376d3d6bff4SMark F. Adams           }
377038e3b61SMark F. Adams         }
378038e3b61SMark F. Adams       }
37911e60469SMark F. Adams       ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr);
38011e60469SMark F. Adams       ierr = VecDestroy(&dest_crd);CHKERRQ(ierr);
381885364a3SMark Adams     }
382a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3830cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
38411e60469SMark F. Adams     /*
3857dae84e0SHong Zhang       Invert for MatCreateSubMatrix
38611e60469SMark F. Adams     */
387a2f3521dSMark F. Adams     ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr);
388a2f3521dSMark F. Adams     ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */
389c5bfad50SMark F. Adams     ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr);
390a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
391a2f3521dSMark F. Adams       ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */
392a2f3521dSMark F. Adams     }
3933cb8563fSToby Isaac     if (Pcolumnperm) {
3943cb8563fSToby Isaac       ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr);
3953cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3963cb8563fSToby Isaac     }
397a2f3521dSMark F. Adams     ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr);
3980cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr);
3990cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
400a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
401a2f3521dSMark F. Adams     {
402a2f3521dSMark F. Adams       Mat       mat;
40390db8557SMark Adams       PetscBool flg;
4047dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr);
40590db8557SMark Adams       ierr = MatGetOption(Cmat, MAT_SPD, &flg);CHKERRQ(ierr);
40690db8557SMark Adams       if (flg) {
40790db8557SMark Adams         ierr = MatSetOption(mat, MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
40890db8557SMark Adams       } else {
40990db8557SMark Adams         ierr = MatGetOption(Cmat, MAT_HERMITIAN, &flg);CHKERRQ(ierr);
41090db8557SMark Adams         if (flg) {
41190db8557SMark Adams           ierr = MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
41290db8557SMark Adams         } else {
41390db8557SMark Adams #if !defined(PETSC_USE_COMPLEX)
41490db8557SMark Adams           ierr = MatGetOption(Cmat, MAT_SYMMETRIC, &flg);CHKERRQ(ierr);
41590db8557SMark Adams           if (flg) {
41690db8557SMark Adams             ierr = MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
41790db8557SMark Adams           }
41890db8557SMark Adams #endif
41990db8557SMark Adams         }
42090db8557SMark Adams       }
421a2f3521dSMark F. Adams       *a_Amat_crs = mat;
422a2f3521dSMark F. Adams     }
423038e3b61SMark F. Adams     ierr = MatDestroy(&Cmat);CHKERRQ(ierr);
424a2f3521dSMark F. Adams 
4250cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr);
42611e60469SMark F. Adams     /* prolongator */
42711e60469SMark F. Adams     {
42811e60469SMark F. Adams       IS       findices;
429a2f3521dSMark F. Adams       PetscInt Istart,Iend;
430a2f3521dSMark F. Adams       Mat      Pnew;
43162294041SBarry Smith 
432a2f3521dSMark F. Adams       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
4330cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
4343b4367a7SBarry Smith       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr);
435c5bfad50SMark F. Adams       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
4367dae84e0SHong Zhang       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
43711e60469SMark F. Adams       ierr = ISDestroy(&findices);CHKERRQ(ierr);
4381a2c6b5cSJunchao Zhang       ierr = MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr);
439c5bfad50SMark F. Adams 
4400cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr);
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;
4541bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: Pinning level %D to the CPU\n",((PetscObject)pc)->prefix,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);
458adf5291fSStefano Zampini       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; 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     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr);
470a2f3521dSMark F. Adams   }
4715b89ad90SMark F. Adams   PetscFunctionReturn(0);
4725b89ad90SMark F. Adams }
4735b89ad90SMark F. Adams 
4744b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
4754b1575e2SStefano Zampini {
4764b1575e2SStefano Zampini   PetscErrorCode ierr;
4774b1575e2SStefano Zampini   const char     *prefix;
4784b1575e2SStefano Zampini   char           addp[32];
4794b1575e2SStefano Zampini   PC_MG          *mg      = (PC_MG*)a_pc->data;
4804b1575e2SStefano Zampini   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4814b1575e2SStefano Zampini 
4824b1575e2SStefano Zampini   PetscFunctionBegin;
4834b1575e2SStefano Zampini   ierr = PCGetOptionsPrefix(a_pc,&prefix);CHKERRQ(ierr);
4841bde8fbeSMark Adams   ierr = PetscInfo(a_pc,"%s: Square Graph on level %D\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1);CHKERRQ(ierr);
4854b1575e2SStefano Zampini   ierr = MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);CHKERRQ(ierr);
4864b1575e2SStefano Zampini   ierr = MatSetOptionsPrefix(*Gmat2,prefix);CHKERRQ(ierr);
4874b1575e2SStefano Zampini   ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);CHKERRQ(ierr);
4884b1575e2SStefano Zampini   ierr = MatAppendOptionsPrefix(*Gmat2,addp);CHKERRQ(ierr);
489b4da3a1bSStefano Zampini   if ((*Gmat2)->structurally_symmetric) {
490b4da3a1bSStefano Zampini     ierr = MatProductSetType(*Gmat2,MATPRODUCT_AB);CHKERRQ(ierr);
491b4da3a1bSStefano Zampini   } else {
4921a2c6b5cSJunchao Zhang     ierr = MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr);
4934b1575e2SStefano Zampini     ierr = MatProductSetType(*Gmat2,MATPRODUCT_AtB);CHKERRQ(ierr);
494b4da3a1bSStefano Zampini   }
4954b1575e2SStefano Zampini   ierr = MatProductSetFromOptions(*Gmat2);CHKERRQ(ierr);
4964555aa8cSStefano Zampini   ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr);
4974b1575e2SStefano Zampini   ierr = MatProductSymbolic(*Gmat2);CHKERRQ(ierr);
4984555aa8cSStefano Zampini   ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);CHKERRQ(ierr);
499b4da3a1bSStefano Zampini   ierr = MatProductClear(*Gmat2);CHKERRQ(ierr);
5004b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
5014b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
5024b1575e2SStefano Zampini   PetscFunctionReturn(0);
5034b1575e2SStefano Zampini }
5044b1575e2SStefano Zampini 
5055b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5065b89ad90SMark F. Adams /*
5075b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5085b89ad90SMark F. Adams                     by setting data structures and options.
5095b89ad90SMark F. Adams 
5105b89ad90SMark F. Adams    Input Parameter:
5115b89ad90SMark F. Adams .  pc - the preconditioner context
5125b89ad90SMark F. Adams 
5135b89ad90SMark F. Adams */
5149d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5155b89ad90SMark F. Adams {
5165b89ad90SMark F. Adams   PetscErrorCode ierr;
5179d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5185b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5192adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
52018c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
5213b4367a7SBarry Smith   MPI_Comm       comm;
522c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
52318c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
52418c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
525a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
526569f4572SMark Adams   MatInfo        info;
527171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
5285ef31b24SMark F. Adams 
5295b89ad90SMark F. Adams   PetscFunctionBegin;
5303b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
531ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
532ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
533dfd5c07aSMark F. Adams 
5348abdc6daSStefano Zampini   if (pc->setupcalled) {
5358abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
536878e152fSMark F. Adams       /* reset everything */
537878e152fSMark F. Adams       ierr = PCReset_MG(pc);CHKERRQ(ierr);
538878e152fSMark F. Adams       pc->setupcalled = 0;
539806fa848SBarry Smith     } else {
54084d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
54103a628feSMark F. Adams       /* just do Galerkin grids */
54258471d46SMark F. Adams       Mat          B,dA,dB;
54358471d46SMark F. Adams 
5449d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
5454555aa8cSStefano Zampini         PetscInt gl;
54658471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
54723ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
54858471d46SMark F. Adams         /* (re)set to get dirty flag */
54923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr);
55058471d46SMark F. Adams 
5514555aa8cSStefano Zampini         for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
5528abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX ;
5538abdc6daSStefano Zampini 
5548abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
55523ee1639SBarry Smith           ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr);
5568abdc6daSStefano Zampini           if (B->product) {
5578abdc6daSStefano Zampini             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
5588abdc6daSStefano Zampini               reuse = MAT_REUSE_MATRIX;
55903a628feSMark F. Adams             }
5608abdc6daSStefano Zampini           }
5618abdc6daSStefano Zampini           if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); }
5628abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
5631bde8fbeSMark Adams             ierr = PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %D\n",((PetscObject)pc)->prefix,level);CHKERRQ(ierr);
5648abdc6daSStefano Zampini           } else {
5651bde8fbeSMark Adams             ierr = PetscInfo(pc,"%s: RAP after first solve, new matrix level %D\n",((PetscObject)pc)->prefix,level);CHKERRQ(ierr);
5668abdc6daSStefano Zampini           }
5674555aa8cSStefano Zampini           ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr);
5688abdc6daSStefano Zampini           ierr = MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);CHKERRQ(ierr);
5694555aa8cSStefano Zampini           ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);CHKERRQ(ierr);
57063b77682SMark Adams           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
57123ee1639SBarry Smith           ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr);
57258471d46SMark F. Adams           dB   = B;
57358471d46SMark F. Adams         }
5745f8cf99dSMark F. Adams       }
575d5280255SMark F. Adams 
576d5280255SMark F. Adams       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
57758471d46SMark F. Adams       PetscFunctionReturn(0);
578eb07cef2SMark F. Adams     }
579878e152fSMark F. Adams   }
580f6536408SMark F. Adams 
581878e152fSMark F. Adams   if (!pc_gamg->data) {
582878e152fSMark F. Adams     if (pc_gamg->orig_data) {
583878e152fSMark F. Adams       ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
5840298fd71SBarry Smith       ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr);
5852fa5cd67SKarl Rupp 
586878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
587878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
588878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5892fa5cd67SKarl Rupp 
590785e854fSJed Brown       ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr);
591878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
592806fa848SBarry Smith     } else {
5932c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5947700e67bSMark Adams       ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr);
5959d5b6da9SMark F. Adams     }
596878e152fSMark F. Adams   }
597878e152fSMark F. Adams 
598878e152fSMark F. Adams   /* cache original data for reuse */
5991c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
600785e854fSJed Brown     ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr);
601878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
602878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
603878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
604878e152fSMark F. Adams   }
605038e3b61SMark F. Adams 
606302f38e8SMark F. Adams   /* get basic dims */
607302f38e8SMark F. Adams   ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr);
608171cca9aSMark Adams   ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr);
60984d3f75bSMark F. Adams 
610569f4572SMark Adams   ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */
611569f4572SMark Adams   nnz0   = info.nz_used;
612569f4572SMark Adams   nnztot = info.nz_used;
6131bde8fbeSMark Adams   ierr = PetscInfo(pc,"%s: level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",((PetscObject)pc)->prefix,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr);
614569f4572SMark Adams 
615a2f3521dSMark F. Adams   /* Get A_i and R_i */
61662294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
6179ab59c8bSMark Adams     pc_gamg->current_level = level;
6182c71b3e2SJacob Faibussowitsch     PetscCheckFalse(level >= PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
6195b89ad90SMark F. Adams     level1 = level + 1;
6200cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
6214555aa8cSStefano Zampini #if defined(GAMG_STAGES)
622a2f3521dSMark F. Adams     ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr);
623b4fbaa2aSMark F. Adams #endif
624c8b0795cSMark F. Adams     { /* construct prolongator */
625725b86d8SJed Brown       Mat              Gmat;
6260cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6277700e67bSMark Adams       Mat              Prol11;
628c8b0795cSMark F. Adams 
6297700e67bSMark Adams       ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr);
6301ab5ffc9SJed Brown       ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr);
6317700e67bSMark Adams       ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr);
632c8b0795cSMark F. Adams 
633a2f3521dSMark F. Adams       /* could have failed to create new level */
634a2f3521dSMark F. Adams       if (Prol11) {
635f7df55f0SStefano Zampini         const char *prefix;
636f7df55f0SStefano Zampini         char       addp[32];
637f7df55f0SStefano Zampini 
6389d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6390298fd71SBarry Smith         ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr);
640a2f3521dSMark F. Adams 
641fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
642c8b0795cSMark F. Adams           /* smooth */
643fd1112cbSBarry Smith           ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr);
644c8b0795cSMark F. Adams         }
645c8b0795cSMark F. Adams 
6460c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
6471b18a24aSMark Adams           PetscInt bs;
6481b18a24aSMark Adams           ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr);
6490a3c815dSMark Adams           ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr);
650ffc955d6SMark F. Adams         }
651ffc955d6SMark F. Adams 
652f7df55f0SStefano Zampini         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
653f7df55f0SStefano Zampini         ierr = MatSetOptionsPrefix(Prol11,prefix);CHKERRQ(ierr);
654c9567895SMark         ierr = PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);CHKERRQ(ierr);
655f7df55f0SStefano Zampini         ierr = MatAppendOptionsPrefix(Prol11,addp);CHKERRQ(ierr);
65691f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
657f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
6581a2c6b5cSJunchao Zhang         ierr = MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);CHKERRQ(ierr);
659f7df55f0SStefano Zampini         ierr = MatSetFromOptions(Prol11);CHKERRQ(ierr);
6604bde40a0SMark Adams         Parr[level1] = Prol11;
6614bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6624bde40a0SMark Adams 
663a2f3521dSMark F. Adams       ierr = MatDestroy(&Gmat);CHKERRQ(ierr);
66441b27cdeSMark F. Adams       ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr);
665a2f3521dSMark F. Adams     } /* construct prolongator scope */
6660cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr);
6677f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
668171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
6691bde8fbeSMark Adams       ierr = PetscInfo(pc,"%s: Stop gridding, level %D\n",((PetscObject)pc)->prefix,level);CHKERRQ(ierr);
6704555aa8cSStefano Zampini #if defined(GAMG_STAGES)
671a90e85d9SMark Adams       ierr = PetscLogStagePop();CHKERRQ(ierr);
672a90e85d9SMark Adams #endif
673c8b0795cSMark F. Adams       break;
674c8b0795cSMark F. Adams     }
6750cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
676171cca9aSMark Adams     ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */
6772c71b3e2SJacob Faibussowitsch     PetscCheckFalse(is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
678171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6790e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
680171cca9aSMark Adams     ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr);
681a2f3521dSMark F. Adams 
6820cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr);
683171cca9aSMark Adams     ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */
684569f4572SMark Adams     ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr);
685569f4572SMark Adams     nnztot += info.nz_used;
6861bde8fbeSMark Adams     ierr = PetscInfo(pc,"%s: %D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",((PetscObject)pc)->prefix,level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);CHKERRQ(ierr);
687569f4572SMark Adams 
6884555aa8cSStefano Zampini #if defined(GAMG_STAGES)
689b4fbaa2aSMark F. Adams     ierr = PetscLogStagePop();CHKERRQ(ierr);
690b4fbaa2aSMark F. Adams #endif
691a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6929ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
6931bde8fbeSMark Adams       ierr =  PetscInfo(pc,"%s: HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",((PetscObject)pc)->prefix,level,M/bs);CHKERRQ(ierr);
694a90e85d9SMark Adams       level++;
695a90e85d9SMark Adams       break;
696a90e85d9SMark Adams     }
697c8b0795cSMark F. Adams   } /* levels */
698c8b0795cSMark F. Adams   ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
699c8b0795cSMark F. Adams 
7001bde8fbeSMark Adams   ierr = PetscInfo(pc,"%s: %D levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0);CHKERRQ(ierr);
7019d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7025b89ad90SMark F. Adams   fine_level       = level;
7030298fd71SBarry Smith   ierr             = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
7045b89ad90SMark F. Adams 
70562294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7060ed2132dSStefano Zampini 
707d5280255SMark F. Adams     /* set default smoothers & set operators */
70862294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
709ffc955d6SMark F. Adams       KSP smoother;
710ffc955d6SMark F. Adams       PC  subpc;
711a2f3521dSMark F. Adams 
7129d5b6da9SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
713f6536408SMark F. Adams       ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
714ffc955d6SMark F. Adams 
715a2f3521dSMark F. Adams       ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
716a2f3521dSMark F. Adams       /* set ops */
71723ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr);
718a2f3521dSMark F. Adams       ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr);
719a2f3521dSMark F. Adams 
720a2f3521dSMark F. Adams       /* set defaults */
7216c9de887SHong Zhang       ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr);
722a2f3521dSMark F. Adams 
7230c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
7240c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
7252d3561bbSSatish Balay         PetscInt sz;
7267a28f3e5SMark Adams         IS       *iss;
727a2f3521dSMark F. Adams 
7282d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7297a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
7300c3bc534SBarry Smith         ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr);
7310a3c815dSMark Adams         ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr);
7320c3bc534SBarry Smith         ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr);
7337f66b68fSMark Adams         if (!sz) {
734ffc955d6SMark F. Adams           IS       is;
7350a3c815dSMark Adams           ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
7367a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr);
737a94c3b12SMark F. Adams           ierr = ISDestroy(&is);CHKERRQ(ierr);
738806fa848SBarry Smith         } else {
739a94c3b12SMark F. Adams           PetscInt kk;
7407a28f3e5SMark Adams           ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr);
741a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
7427a28f3e5SMark Adams             ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr);
743a94c3b12SMark F. Adams           }
7447a28f3e5SMark Adams           ierr = PetscFree(iss);CHKERRQ(ierr);
745ffc955d6SMark F. Adams         }
7460298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
747ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
748806fa848SBarry Smith       } else {
7495f7df010SMark Adams         ierr = PCSetType(subpc, PCJACOBI);CHKERRQ(ierr);
750ffc955d6SMark F. Adams       }
751d5280255SMark F. Adams     }
752d5280255SMark F. Adams     {
753d5280255SMark F. Adams       /* coarse grid */
754d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
755d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
7560ed2132dSStefano Zampini 
757d5280255SMark F. Adams       ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
75823ee1639SBarry Smith       ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr);
759cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
760d5280255SMark F. Adams         ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);
761d5280255SMark F. Adams         ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr);
762d5280255SMark F. Adams         ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr);
763d5280255SMark F. Adams         ierr = PCSetUp(subpc);CHKERRQ(ierr);
76471959b99SBarry Smith         ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr);
7652c71b3e2SJacob Faibussowitsch         PetscCheckFalse(ii != 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
766d5280255SMark F. Adams         ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr);
767d5280255SMark F. Adams         ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr);
7689dbfc187SHong Zhang         ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
7692fb0b348SMark F. Adams         ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
77008e36f19SMark Adams         ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr);
771d5280255SMark F. Adams       }
772cf8ae1d3SMark Adams     }
773d5280255SMark F. Adams 
774d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
775d5280255SMark F. Adams     ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
776e55864a3SBarry Smith     ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr);
777d5280255SMark F. Adams     ierr = PetscOptionsEnd();CHKERRQ(ierr);
77869aca0b8SBarry Smith     ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
779d5280255SMark F. Adams 
78018c3aa7eSMark     /* setup cheby eigen estimates from SA */
7817e6512fdSJed Brown     if (pc_gamg->use_sa_esteig) {
78218c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
78318c3aa7eSMark         KSP       smoother;
78418c3aa7eSMark         PetscBool ischeb;
7850ed2132dSStefano Zampini 
78618c3aa7eSMark         ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr);
78718c3aa7eSMark         ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr);
78818c3aa7eSMark         if (ischeb) {
78918c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
7900ed2132dSStefano Zampini 
7912de708cbSJed Brown           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
792efe053fcSJed Brown           if (mg->max_eigen_DinvA[level] > 0) {
7937e6512fdSJed Brown             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
7947e6512fdSJed Brown             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
79518c3aa7eSMark             PetscReal emax,emin;
7960ed2132dSStefano Zampini 
79718c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
79818c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
7991bde8fbeSMark Adams             ierr = PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin);CHKERRQ(ierr);
8000a94a983SJed Brown             cheb->emin_provided = emin;
8010a94a983SJed Brown             cheb->emax_provided = emax;
80218c3aa7eSMark           }
80318c3aa7eSMark         }
80418c3aa7eSMark       }
80518c3aa7eSMark     }
8060ed2132dSStefano Zampini 
8070ed2132dSStefano Zampini     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8080ed2132dSStefano Zampini 
809d5280255SMark F. Adams     /* clean up */
810d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
811587fa25dSMark F. Adams       ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr);
812587fa25dSMark F. Adams       ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr);
8135b89ad90SMark F. Adams     }
814806fa848SBarry Smith   } else {
8155f8cf99dSMark F. Adams     KSP smoother;
8160ed2132dSStefano Zampini 
8171bde8fbeSMark Adams     ierr = PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr);
8189d5b6da9SMark F. Adams     ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr);
81923ee1639SBarry Smith     ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr);
8205f8cf99dSMark F. Adams     ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr);
8219d5b6da9SMark F. Adams     ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
8225f8cf99dSMark F. Adams   }
8235b89ad90SMark F. Adams   PetscFunctionReturn(0);
8245b89ad90SMark F. Adams }
8255b89ad90SMark F. Adams 
826eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8275b89ad90SMark F. Adams /*
8285b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8295b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8305b89ad90SMark F. Adams 
8315b89ad90SMark F. Adams    Input Parameter:
8325b89ad90SMark F. Adams .  pc - the preconditioner context
8335b89ad90SMark F. Adams 
8345b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8355b89ad90SMark F. Adams */
8365b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8375b89ad90SMark F. Adams {
8385b89ad90SMark F. Adams   PetscErrorCode ierr;
8395b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8405b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8415b89ad90SMark F. Adams 
8425b89ad90SMark F. Adams   PetscFunctionBegin;
8435b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
8449b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8459b8ffb57SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
8469b8ffb57SJed Brown   }
8471ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr);
8481ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
8495b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
8505b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8515b89ad90SMark F. Adams   PetscFunctionReturn(0);
8525b89ad90SMark F. Adams }
8535b89ad90SMark F. Adams 
854676e1743SMark F. Adams /*@
855cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
856676e1743SMark F. Adams 
8571cc46a46SBarry Smith    Logically Collective on PC
858676e1743SMark F. Adams 
859676e1743SMark F. Adams    Input Parameters:
8601cc46a46SBarry Smith +  pc - the preconditioner context
8611cc46a46SBarry Smith -  n - the number of equations
862676e1743SMark F. Adams 
863676e1743SMark F. Adams    Options Database Key:
864*147403d9SBarry Smith .  -pc_gamg_process_eq_limit <limit> - set the limit
865676e1743SMark F. Adams 
86695452b02SPatrick Sanan    Notes:
86795452b02SPatrick 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
868cab9ed1eSBarry Smith     that has degrees of freedom
869cab9ed1eSBarry Smith 
870676e1743SMark F. Adams    Level: intermediate
871676e1743SMark F. Adams 
872c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
873676e1743SMark F. Adams @*/
874676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
875676e1743SMark F. Adams {
876676e1743SMark F. Adams   PetscErrorCode ierr;
877676e1743SMark F. Adams 
878676e1743SMark F. Adams   PetscFunctionBegin;
879676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
880676e1743SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
881676e1743SMark F. Adams   PetscFunctionReturn(0);
882676e1743SMark F. Adams }
883676e1743SMark F. Adams 
8841e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
885676e1743SMark F. Adams {
886c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
887c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
888676e1743SMark F. Adams 
889676e1743SMark F. Adams   PetscFunctionBegin;
8909d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
891676e1743SMark F. Adams   PetscFunctionReturn(0);
892676e1743SMark F. Adams }
893676e1743SMark F. Adams 
894389730f3SMark F. Adams /*@
895cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
896389730f3SMark F. Adams 
897389730f3SMark F. Adams  Collective on PC
898389730f3SMark F. Adams 
899389730f3SMark F. Adams    Input Parameters:
9001cc46a46SBarry Smith +  pc - the preconditioner context
9011cc46a46SBarry Smith -  n - maximum number of equations to aim for
902389730f3SMark F. Adams 
903389730f3SMark F. Adams    Options Database Key:
904*147403d9SBarry Smith .  -pc_gamg_coarse_eq_limit <limit> - set the limit
905389730f3SMark F. Adams 
90674329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
90774329af1SBarry Smith      has less than 1000 unknowns.
90874329af1SBarry Smith 
909389730f3SMark F. Adams    Level: intermediate
910389730f3SMark F. Adams 
911c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
912389730f3SMark F. Adams @*/
913389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
914389730f3SMark F. Adams {
915389730f3SMark F. Adams   PetscErrorCode ierr;
916389730f3SMark F. Adams 
917389730f3SMark F. Adams   PetscFunctionBegin;
918389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
919389730f3SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
920389730f3SMark F. Adams   PetscFunctionReturn(0);
921389730f3SMark F. Adams }
922389730f3SMark F. Adams 
9231e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
924389730f3SMark F. Adams {
925389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
926389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
927389730f3SMark F. Adams 
928389730f3SMark F. Adams   PetscFunctionBegin;
9299d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
930389730f3SMark F. Adams   PetscFunctionReturn(0);
931389730f3SMark F. Adams }
932389730f3SMark F. Adams 
933676e1743SMark F. Adams /*@
934cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
935676e1743SMark F. Adams 
936676e1743SMark F. Adams    Collective on PC
937676e1743SMark F. Adams 
938676e1743SMark F. Adams    Input Parameters:
9391cc46a46SBarry Smith +  pc - the preconditioner context
9401cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
941676e1743SMark F. Adams 
942676e1743SMark F. Adams    Options Database Key:
943*147403d9SBarry Smith .  -pc_gamg_repartition <true,false> - turn on the repartitioning
944676e1743SMark F. Adams 
94595452b02SPatrick Sanan    Notes:
94695452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
947cab9ed1eSBarry Smith 
948676e1743SMark F. Adams    Level: intermediate
949676e1743SMark F. Adams 
950676e1743SMark F. Adams @*/
951cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
952676e1743SMark F. Adams {
953676e1743SMark F. Adams   PetscErrorCode ierr;
954676e1743SMark F. Adams 
955676e1743SMark F. Adams   PetscFunctionBegin;
956676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
957cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
958676e1743SMark F. Adams   PetscFunctionReturn(0);
959676e1743SMark F. Adams }
960676e1743SMark F. Adams 
961cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
962676e1743SMark F. Adams {
963c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
964c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
965676e1743SMark F. Adams 
966676e1743SMark F. Adams   PetscFunctionBegin;
9679d5b6da9SMark F. Adams   pc_gamg->repart = n;
968676e1743SMark F. Adams   PetscFunctionReturn(0);
969676e1743SMark F. Adams }
970676e1743SMark F. Adams 
971dfd5c07aSMark F. Adams /*@
9727e6512fdSJed Brown    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother
97318c3aa7eSMark 
97418c3aa7eSMark    Collective on PC
97518c3aa7eSMark 
97618c3aa7eSMark    Input Parameters:
97718c3aa7eSMark +  pc - the preconditioner context
97818c3aa7eSMark -  n - number of its
97918c3aa7eSMark 
98018c3aa7eSMark    Options Database Key:
981*147403d9SBarry Smith .  -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
98218c3aa7eSMark 
98318c3aa7eSMark    Notes:
9847e6512fdSJed Brown    Smoothed aggregation constructs the smoothed prolongator $P = (I - \omega D^{-1} A) T$ where $T$ is the tentative prolongator and $D$ is the diagonal of $A$.
9857e6512fdSJed Brown    Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
9867e6512fdSJed Brown    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused.
987efe053fcSJed Brown    This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used.
988efe053fcSJed Brown    It became default in PETSc 3.17.
98918c3aa7eSMark 
9907e6512fdSJed Brown    Level: advanced
99118c3aa7eSMark 
99273f7197eSJed Brown .seealso: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet()
99318c3aa7eSMark @*/
99418c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
99518c3aa7eSMark {
99618c3aa7eSMark   PetscErrorCode ierr;
99718c3aa7eSMark 
99818c3aa7eSMark   PetscFunctionBegin;
99918c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
100018c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
100118c3aa7eSMark   PetscFunctionReturn(0);
100218c3aa7eSMark }
100318c3aa7eSMark 
10040ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
100518c3aa7eSMark {
100618c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
100718c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
100818c3aa7eSMark 
100918c3aa7eSMark   PetscFunctionBegin;
10107e6512fdSJed Brown   pc_gamg->use_sa_esteig = n;
101118c3aa7eSMark   PetscFunctionReturn(0);
101218c3aa7eSMark }
101318c3aa7eSMark 
101418c3aa7eSMark /*@
101518c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
101618c3aa7eSMark 
101718c3aa7eSMark    Collective on PC
101818c3aa7eSMark 
101918c3aa7eSMark    Input Parameters:
102018c3aa7eSMark +  pc - the preconditioner context
102118c3aa7eSMark -  emax - max eigenvalue
102218c3aa7eSMark .  emin - min eigenvalue
102318c3aa7eSMark 
102418c3aa7eSMark    Options Database Key:
1025*147403d9SBarry Smith .  -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
102618c3aa7eSMark 
102718c3aa7eSMark    Level: intermediate
102818c3aa7eSMark 
102973f7197eSJed Brown .seealso: PCGAMGSetUseSAEstEig()
103018c3aa7eSMark @*/
103118c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
103218c3aa7eSMark {
103318c3aa7eSMark   PetscErrorCode ierr;
103418c3aa7eSMark 
103518c3aa7eSMark   PetscFunctionBegin;
103618c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
103718c3aa7eSMark   ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr);
103818c3aa7eSMark   PetscFunctionReturn(0);
103918c3aa7eSMark }
104041ffd417SStefano Zampini 
104118c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
104218c3aa7eSMark {
104318c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
104418c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
104518c3aa7eSMark 
104618c3aa7eSMark   PetscFunctionBegin;
10472c71b3e2SJacob Faibussowitsch   PetscCheckFalse(emax <= emin,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
10482c71b3e2SJacob Faibussowitsch   PetscCheckFalse(emax*emin <= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
104918c3aa7eSMark   pc_gamg->emax = emax;
105018c3aa7eSMark   pc_gamg->emin = emin;
105118c3aa7eSMark 
105218c3aa7eSMark   PetscFunctionReturn(0);
105318c3aa7eSMark }
105418c3aa7eSMark 
105518c3aa7eSMark /*@
1056cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1057dfd5c07aSMark F. Adams 
1058dfd5c07aSMark F. Adams    Collective on PC
1059dfd5c07aSMark F. Adams 
1060dfd5c07aSMark F. Adams    Input Parameters:
10611cc46a46SBarry Smith +  pc - the preconditioner context
10621cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1063dfd5c07aSMark F. Adams 
1064dfd5c07aSMark F. Adams    Options Database Key:
1065*147403d9SBarry Smith .  -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1066dfd5c07aSMark F. Adams 
1067dfd5c07aSMark F. Adams    Level: intermediate
1068dfd5c07aSMark F. Adams 
106995452b02SPatrick Sanan    Notes:
1070*147403d9SBarry Smith     May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1071cab9ed1eSBarry Smith     rebuilding the preconditioner quicker.
1072cab9ed1eSBarry Smith 
1073dfd5c07aSMark F. Adams @*/
10741cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1075dfd5c07aSMark F. Adams {
1076dfd5c07aSMark F. Adams   PetscErrorCode ierr;
1077dfd5c07aSMark F. Adams 
1078dfd5c07aSMark F. Adams   PetscFunctionBegin;
1079dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10801cc46a46SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
1081dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1082dfd5c07aSMark F. Adams }
1083dfd5c07aSMark F. Adams 
10841cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1085dfd5c07aSMark F. Adams {
1086dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1087dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1088dfd5c07aSMark F. Adams 
1089dfd5c07aSMark F. Adams   PetscFunctionBegin;
1090dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1091dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1092dfd5c07aSMark F. Adams }
1093dfd5c07aSMark F. Adams 
1094ffc955d6SMark F. Adams /*@
1095cab9ed1eSBarry 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.
1096ffc955d6SMark F. Adams 
1097ffc955d6SMark F. Adams    Collective on PC
1098ffc955d6SMark F. Adams 
1099ffc955d6SMark F. Adams    Input Parameters:
1100cab9ed1eSBarry Smith +  pc - the preconditioner context
1101cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1102ffc955d6SMark F. Adams 
1103ffc955d6SMark F. Adams    Options Database Key:
1104*147403d9SBarry Smith .  -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1105ffc955d6SMark F. Adams 
1106ffc955d6SMark F. Adams    Level: intermediate
1107ffc955d6SMark F. Adams 
1108ffc955d6SMark F. Adams @*/
1109cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1110ffc955d6SMark F. Adams {
1111ffc955d6SMark F. Adams   PetscErrorCode ierr;
1112ffc955d6SMark F. Adams 
1113ffc955d6SMark F. Adams   PetscFunctionBegin;
1114ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1115cab9ed1eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1116ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1117ffc955d6SMark F. Adams }
1118ffc955d6SMark F. Adams 
1119cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1120ffc955d6SMark F. Adams {
1121ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1122ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1123ffc955d6SMark F. Adams 
1124ffc955d6SMark F. Adams   PetscFunctionBegin;
1125cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1126ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1127ffc955d6SMark F. Adams }
1128ffc955d6SMark F. Adams 
1129171cca9aSMark Adams /*@
1130cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1131171cca9aSMark Adams 
1132171cca9aSMark Adams    Collective on PC
1133171cca9aSMark Adams 
1134171cca9aSMark Adams    Input Parameters:
1135171cca9aSMark Adams +  pc - the preconditioner context
1136cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1137171cca9aSMark Adams 
1138171cca9aSMark Adams    Options Database Key:
1139*147403d9SBarry Smith .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1140171cca9aSMark Adams 
1141171cca9aSMark Adams    Level: intermediate
1142171cca9aSMark Adams 
114339d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1144171cca9aSMark Adams @*/
1145171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1146171cca9aSMark Adams {
1147171cca9aSMark Adams   PetscErrorCode ierr;
1148171cca9aSMark Adams 
1149171cca9aSMark Adams   PetscFunctionBegin;
1150171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1151171cca9aSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1152171cca9aSMark Adams   PetscFunctionReturn(0);
1153171cca9aSMark Adams }
1154171cca9aSMark Adams 
1155171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1156171cca9aSMark Adams {
1157171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1158171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1159171cca9aSMark Adams 
1160171cca9aSMark Adams   PetscFunctionBegin;
1161171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1162ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1163ffc955d6SMark F. Adams }
1164ffc955d6SMark F. Adams 
11654ef23d27SMark F. Adams /*@
1166ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1167ce7c7f2fSMark Adams 
1168ce7c7f2fSMark Adams    Collective on PC
1169ce7c7f2fSMark Adams 
1170ce7c7f2fSMark Adams    Input Parameters:
1171ce7c7f2fSMark Adams +  pc - the preconditioner context
1172ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1173ce7c7f2fSMark Adams 
1174ce7c7f2fSMark Adams    Options Database Key:
1175*147403d9SBarry Smith .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1176ce7c7f2fSMark Adams 
1177ce7c7f2fSMark Adams    Level: intermediate
1178ce7c7f2fSMark Adams 
117939d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1180ce7c7f2fSMark Adams @*/
1181ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1182ce7c7f2fSMark Adams {
1183ce7c7f2fSMark Adams   PetscErrorCode ierr;
1184ce7c7f2fSMark Adams 
1185ce7c7f2fSMark Adams   PetscFunctionBegin;
1186ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1187ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1188ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1189ce7c7f2fSMark Adams }
1190ce7c7f2fSMark Adams 
1191ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1192ce7c7f2fSMark Adams {
1193ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1194ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1195ce7c7f2fSMark Adams 
1196ce7c7f2fSMark Adams   PetscFunctionBegin;
1197ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1198ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1199ce7c7f2fSMark Adams }
1200ce7c7f2fSMark Adams 
1201ce7c7f2fSMark Adams /*@
1202*147403d9SBarry Smith    PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1203ce7c7f2fSMark Adams 
1204ce7c7f2fSMark Adams    Collective on PC
1205ce7c7f2fSMark Adams 
1206ce7c7f2fSMark Adams    Input Parameters:
1207ce7c7f2fSMark Adams +  pc - the preconditioner context
1208ce7c7f2fSMark Adams -  flg - Layout type
1209ce7c7f2fSMark Adams 
1210ce7c7f2fSMark Adams    Options Database Key:
1211*147403d9SBarry Smith .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1212ce7c7f2fSMark Adams 
1213ce7c7f2fSMark Adams    Level: intermediate
1214ce7c7f2fSMark Adams 
121539d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1216ce7c7f2fSMark Adams @*/
1217ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1218ce7c7f2fSMark Adams {
1219ce7c7f2fSMark Adams   PetscErrorCode ierr;
1220ce7c7f2fSMark Adams 
1221ce7c7f2fSMark Adams   PetscFunctionBegin;
1222ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1223ce7c7f2fSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr);
1224ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1225ce7c7f2fSMark Adams }
1226ce7c7f2fSMark Adams 
1227ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1228ce7c7f2fSMark Adams {
1229ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1230ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1231ce7c7f2fSMark Adams 
1232ce7c7f2fSMark Adams   PetscFunctionBegin;
1233ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1234ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1235ce7c7f2fSMark Adams }
1236ce7c7f2fSMark Adams 
1237ce7c7f2fSMark Adams /*@
12381cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
12394ef23d27SMark F. Adams 
12404ef23d27SMark F. Adams    Not collective on PC
12414ef23d27SMark F. Adams 
12424ef23d27SMark F. Adams    Input Parameters:
12431cc46a46SBarry Smith +  pc - the preconditioner
12441cc46a46SBarry Smith -  n - the maximum number of levels to use
12454ef23d27SMark F. Adams 
12464ef23d27SMark F. Adams    Options Database Key:
1247*147403d9SBarry Smith .  -pc_mg_levels <n> - set the maximum number of levels to allow
12484ef23d27SMark F. Adams 
12494ef23d27SMark F. Adams    Level: intermediate
12504ef23d27SMark F. Adams 
12514ef23d27SMark F. Adams @*/
12524ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12534ef23d27SMark F. Adams {
12544ef23d27SMark F. Adams   PetscErrorCode ierr;
12554ef23d27SMark F. Adams 
12564ef23d27SMark F. Adams   PetscFunctionBegin;
12574ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12584ef23d27SMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
12594ef23d27SMark F. Adams   PetscFunctionReturn(0);
12604ef23d27SMark F. Adams }
12614ef23d27SMark F. Adams 
12621e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12634ef23d27SMark F. Adams {
12644ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12654ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12664ef23d27SMark F. Adams 
12674ef23d27SMark F. Adams   PetscFunctionBegin;
12689d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12694ef23d27SMark F. Adams   PetscFunctionReturn(0);
12704ef23d27SMark F. Adams }
12714ef23d27SMark F. Adams 
12723542efc5SMark F. Adams /*@
12733542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12743542efc5SMark F. Adams 
12753542efc5SMark F. Adams    Not collective on PC
12763542efc5SMark F. Adams 
12773542efc5SMark F. Adams    Input Parameters:
12781cc46a46SBarry Smith +  pc - the preconditioner context
1279c9567895SMark .  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
1280055c8bd0SJed Brown -  n - number of threshold values provided in array
12813542efc5SMark F. Adams 
12823542efc5SMark F. Adams    Options Database Key:
1283*147403d9SBarry Smith .  -pc_gamg_threshold <threshold> - the threshold to drop edges
12843542efc5SMark F. Adams 
128595452b02SPatrick Sanan    Notes:
1286af3c827dSMark 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.
1287af3c827dSMark 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.
1288cab9ed1eSBarry Smith 
1289055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1290055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1291055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1292055c8bd0SJed Brown 
12933542efc5SMark F. Adams    Level: intermediate
12943542efc5SMark F. Adams 
1295af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
12963542efc5SMark F. Adams @*/
1297c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
12983542efc5SMark F. Adams {
12993542efc5SMark F. Adams   PetscErrorCode ierr;
13003542efc5SMark F. Adams 
13013542efc5SMark F. Adams   PetscFunctionBegin;
13023542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1303055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
1304c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr);
13053542efc5SMark F. Adams   PetscFunctionReturn(0);
13063542efc5SMark F. Adams }
13073542efc5SMark F. Adams 
1308c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
13093542efc5SMark F. Adams {
1310c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1311c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1312c1eae691SMark Adams   PetscInt i;
1313c1eae691SMark Adams   PetscFunctionBegin;
1314055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1315055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1316c1eae691SMark Adams   PetscFunctionReturn(0);
1317c1eae691SMark Adams }
1318c1eae691SMark Adams 
1319c1eae691SMark Adams /*@
1320*147403d9SBarry Smith    PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids
1321c9567895SMark 
1322c9567895SMark    Collective on PC
1323c9567895SMark 
1324c9567895SMark    Input Parameters:
1325c9567895SMark +  pc - the preconditioner context
1326c9567895SMark .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1327c9567895SMark -  n - number of values provided in array
1328c9567895SMark 
1329c9567895SMark    Options Database Key:
1330*147403d9SBarry Smith .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1331c9567895SMark 
1332c9567895SMark    Level: intermediate
1333c9567895SMark 
1334c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1335c9567895SMark @*/
1336c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1337c9567895SMark {
1338c9567895SMark   PetscErrorCode ierr;
1339c9567895SMark 
1340c9567895SMark   PetscFunctionBegin;
1341c9567895SMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1342c9567895SMark   if (n) PetscValidIntPointer(v,2);
1343c9567895SMark   ierr = PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));CHKERRQ(ierr);
1344c9567895SMark   PetscFunctionReturn(0);
1345c9567895SMark }
1346c9567895SMark 
1347c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1348c9567895SMark {
1349c9567895SMark   PC_MG   *mg      = (PC_MG*)pc->data;
1350c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1351c9567895SMark   PetscInt i;
1352c9567895SMark   PetscFunctionBegin;
1353c9567895SMark   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1354c9567895SMark   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1355c9567895SMark   PetscFunctionReturn(0);
1356c9567895SMark }
1357c9567895SMark 
1358c9567895SMark /*@
1359c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1360c1eae691SMark Adams 
1361c1eae691SMark Adams    Not collective on PC
1362c1eae691SMark Adams 
1363c1eae691SMark Adams    Input Parameters:
1364c1eae691SMark Adams +  pc - the preconditioner context
1365*147403d9SBarry Smith -  scale - the threshold value reduction, usually < 1.0
1366c1eae691SMark Adams 
1367c1eae691SMark Adams    Options Database Key:
1368*147403d9SBarry Smith .  -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1369c1eae691SMark Adams 
1370055c8bd0SJed Brown    Notes:
1371055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1372055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1373055c8bd0SJed Brown 
1374c1eae691SMark Adams    Level: advanced
1375c1eae691SMark Adams 
1376055c8bd0SJed Brown .seealso: PCGAMGSetThreshold()
1377c1eae691SMark Adams @*/
1378c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1379c1eae691SMark Adams {
1380c1eae691SMark Adams   PetscErrorCode ierr;
13813542efc5SMark F. Adams 
13823542efc5SMark F. Adams   PetscFunctionBegin;
1383c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1384c1eae691SMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr);
1385c1eae691SMark Adams   PetscFunctionReturn(0);
1386c1eae691SMark Adams }
1387c1eae691SMark Adams 
1388c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1389c1eae691SMark Adams {
1390c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1391c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1392c1eae691SMark Adams   PetscFunctionBegin;
1393c1eae691SMark Adams   pc_gamg->threshold_scale = v;
13943542efc5SMark F. Adams   PetscFunctionReturn(0);
13953542efc5SMark F. Adams }
13963542efc5SMark F. Adams 
1397e20c40e8SBarry Smith /*@C
1398c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1399676e1743SMark F. Adams 
1400676e1743SMark F. Adams    Collective on PC
1401676e1743SMark F. Adams 
1402676e1743SMark F. Adams    Input Parameters:
1403c60c7ad4SBarry Smith +  pc - the preconditioner context
1404c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1405676e1743SMark F. Adams 
1406676e1743SMark F. Adams    Options Database Key:
1407cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1408676e1743SMark F. Adams 
1409676e1743SMark F. Adams    Level: intermediate
1410676e1743SMark F. Adams 
1411cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1412676e1743SMark F. Adams @*/
141319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1414676e1743SMark F. Adams {
1415676e1743SMark F. Adams   PetscErrorCode ierr;
1416676e1743SMark F. Adams 
1417676e1743SMark F. Adams   PetscFunctionBegin;
1418676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1419806fa848SBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
1420676e1743SMark F. Adams   PetscFunctionReturn(0);
1421676e1743SMark F. Adams }
1422676e1743SMark F. Adams 
1423e20c40e8SBarry Smith /*@C
1424c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1425c60c7ad4SBarry Smith 
1426c60c7ad4SBarry Smith    Collective on PC
1427c60c7ad4SBarry Smith 
1428c60c7ad4SBarry Smith    Input Parameter:
1429c60c7ad4SBarry Smith .  pc - the preconditioner context
1430c60c7ad4SBarry Smith 
1431c60c7ad4SBarry Smith    Output Parameter:
1432c60c7ad4SBarry Smith .  type - the type of algorithm used
1433c60c7ad4SBarry Smith 
1434c60c7ad4SBarry Smith    Level: intermediate
1435c60c7ad4SBarry Smith 
14361c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1437c60c7ad4SBarry Smith @*/
1438c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1439c60c7ad4SBarry Smith {
1440c60c7ad4SBarry Smith   PetscErrorCode ierr;
1441c60c7ad4SBarry Smith 
1442c60c7ad4SBarry Smith   PetscFunctionBegin;
1443c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1444c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
1445c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1446c60c7ad4SBarry Smith }
1447c60c7ad4SBarry Smith 
1448c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1449c60c7ad4SBarry Smith {
1450c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1451c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1452c60c7ad4SBarry Smith 
1453c60c7ad4SBarry Smith   PetscFunctionBegin;
1454c60c7ad4SBarry Smith   *type = pc_gamg->type;
1455c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1456c60c7ad4SBarry Smith }
1457c60c7ad4SBarry Smith 
14581e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1459676e1743SMark F. Adams {
14609d5b6da9SMark F. Adams   PetscErrorCode ierr,(*r)(PC);
14611ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
14621ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1463676e1743SMark F. Adams 
1464676e1743SMark F. Adams   PetscFunctionBegin;
1465c60c7ad4SBarry Smith   pc_gamg->type = type;
14661c9cd337SJed Brown   ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr);
14672c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
14681ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
14691ab5ffc9SJed Brown     ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr);
14701ab5ffc9SJed Brown     ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr);
1471e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14723ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
14733ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
14743ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
14753ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
14763ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
14773ae0bb68SMark Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
14783ae0bb68SMark Adams     pc_gamg->data_sz = 0;
14791ab5ffc9SJed Brown   }
14801ab5ffc9SJed Brown   ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr);
14811ab5ffc9SJed Brown   ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr);
14829d5b6da9SMark F. Adams   ierr = (*r)(pc);CHKERRQ(ierr);
1483676e1743SMark F. Adams   PetscFunctionReturn(0);
1484676e1743SMark F. Adams }
1485676e1743SMark F. Adams 
14865adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
14875adeb434SBarry Smith {
1488c1eae691SMark Adams   PetscErrorCode ierr,i;
14895adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
14905adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1491e7d4b4cbSMark Adams   PetscReal       gc=0, oc=0;
149290db8557SMark Adams 
14935adeb434SBarry Smith   PetscFunctionBegin;
14945adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");CHKERRQ(ierr);
1495459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");CHKERRQ(ierr);
1496b3e187dcSStefano Zampini   for (i=0;i<mg->nlevels; i++) {
1497c1eae691SMark Adams     ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr);
1498c1eae691SMark Adams   }
1499459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1500459726d8SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr);
1501cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
1502cab9ed1eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr);
1503cab9ed1eSBarry Smith   }
1504171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
1505171cca9aSMark Adams     ierr = PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr);
1506171cca9aSMark Adams   }
1507ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1508ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
1509ce7c7f2fSMark Adams     /* ierr = PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */
1510ce7c7f2fSMark Adams   }
1511ce7c7f2fSMark Adams #endif
1512ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1513ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */
1514ce7c7f2fSMark Adams   /* } else { */
1515ce7c7f2fSMark Adams   /*   ierr = PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */
1516ce7c7f2fSMark Adams   /* } */
15175adeb434SBarry Smith   if (pc_gamg->ops->view) {
15185adeb434SBarry Smith     ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr);
15195adeb434SBarry Smith   }
1520e7d4b4cbSMark Adams   ierr = PCMGGetGridComplexity(pc,&gc,&oc);CHKERRQ(ierr);
1521e7d4b4cbSMark Adams   ierr = PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g    operator = %g\n",gc,oc);CHKERRQ(ierr);
15225adeb434SBarry Smith   PetscFunctionReturn(0);
15235adeb434SBarry Smith }
15245adeb434SBarry Smith 
15254416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
15265b89ad90SMark F. Adams {
1527676e1743SMark F. Adams   PetscErrorCode ierr;
1528676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1529676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
15307e6512fdSJed Brown   PetscBool      flag;
15313b4367a7SBarry Smith   MPI_Comm       comm;
153218c3aa7eSMark   char           prefix[256],tname[32];
1533c1eae691SMark Adams   PetscInt       i,n;
153414a9496bSBarry Smith   const char     *pcpre;
15350a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
15365b89ad90SMark F. Adams   PetscFunctionBegin;
15373b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1538e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr);
15391a1c1e04SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr);
1540bd94a7aaSJed Brown   if (flag) {
1541bd94a7aaSJed Brown     ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr);
15421ab5ffc9SJed Brown   }
1543cab9ed1eSBarry Smith   ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr);
1544efe053fcSJed Brown   ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",pc_gamg->use_sa_esteig,&pc_gamg->use_sa_esteig,NULL);CHKERRQ(ierr);
15451cc46a46SBarry Smith   ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr);
1546a303c832SJed 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);
1547cf8ae1d3SMark 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);
1548ce7c7f2fSMark 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);
1549a0095786SMark   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);
155094ae4db5SBarry 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);
155194ae4db5SBarry 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);
1552a303c832SJed 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);
155318c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
1554c1eae691SMark Adams   ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr);
155518c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1556efd3c5ceSMark Adams     if (!flag) n = 1;
1557c1eae691SMark Adams     i = n;
155818c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1559c1eae691SMark Adams   }
1560c9567895SMark   n = PETSC_MG_MAXLEVELS;
1561c9567895SMark   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);
1562c9567895SMark   if (!flag) i = 0;
1563c9567895SMark   else i = n;
1564c9567895SMark   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
156594ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr);
156618c3aa7eSMark   {
156718c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
156818c3aa7eSMark     n = 2;
156918c3aa7eSMark     ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr);
157018c3aa7eSMark     if (flag) {
15712c71b3e2SJacob Faibussowitsch       PetscCheckFalse(n != 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
157218c3aa7eSMark       ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr);
157318c3aa7eSMark     }
157418c3aa7eSMark   }
1575b7cbab4eSMark Adams   /* set options for subtype */
1576e55864a3SBarry Smith   if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);}
157718c3aa7eSMark 
157814a9496bSBarry Smith   ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr);
157914a9496bSBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr);
1580676e1743SMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
15815b89ad90SMark F. Adams   PetscFunctionReturn(0);
15825b89ad90SMark F. Adams }
15835b89ad90SMark F. Adams 
15845b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
15855b89ad90SMark F. Adams /*MC
15861cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
15875b89ad90SMark F. Adams 
1588280d9858SJed Brown    Options Database Keys:
1589cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1590cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1591cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1592cab9ed1eSBarry 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
1593cab9ed1eSBarry 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>
1594cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1595cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
15966008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1597c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1598cab9ed1eSBarry Smith 
1599cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1600cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1601cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1602cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1603cab9ed1eSBarry Smith 
1604db9745e2SBarry Smith    Multigrid options:
1605db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1606db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1607db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1608cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
16095b89ad90SMark F. Adams 
161095452b02SPatrick Sanan   Notes:
161195452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1612db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1613db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1614db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
16151cc46a46SBarry Smith 
16165b89ad90SMark F. Adams   Level: intermediate
1617280d9858SJed Brown 
16181cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
161973f7197eSJed Brown            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig()
16205b89ad90SMark F. Adams M*/
1621b2573a8aSBarry Smith 
16228cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
16235b89ad90SMark F. Adams {
1624c1eae691SMark Adams   PetscErrorCode ierr,i;
16255b89ad90SMark F. Adams   PC_GAMG        *pc_gamg;
16265b89ad90SMark F. Adams   PC_MG          *mg;
16275b89ad90SMark F. Adams 
16285b89ad90SMark F. Adams   PetscFunctionBegin;
16291c1aac46SBarry Smith    /* register AMG type */
16301c1aac46SBarry Smith   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
16311c1aac46SBarry Smith 
16325b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
16331c1aac46SBarry Smith   ierr = PCSetType(pc, PCMG);CHKERRQ(ierr);
16345b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr);
16355b89ad90SMark F. Adams 
16365b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
1637b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr);
163869aca0b8SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
16395b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
16405b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
16415b89ad90SMark F. Adams 
1642b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr);
16431ab5ffc9SJed Brown 
16449d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
16459d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
16460a545947SLisandro Dalcin   pc_gamg->data    = NULL;
16475b89ad90SMark F. Adams 
16489d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
16495b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
16505b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
16515b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
16525b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
16535adeb434SBarry Smith   mg->view                = PCView_GAMG;
16545b89ad90SMark F. Adams 
165597d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
165697d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1657bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr);
1658bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr);
1659cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr);
166018c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr);
166118c3aa7eSMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr);
16621cc46a46SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr);
1663cab9ed1eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr);
1664171cca9aSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr);
1665ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr);
1666ce7c7f2fSMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr);
1667bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr);
1668c9567895SMark   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);CHKERRQ(ierr);
1669c1eae691SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr);
1670bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr);
1671c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr);
1672bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr);
16739d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1674d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
16750c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1676171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1677a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1678a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1679038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
168025a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
168118c3aa7eSMark   for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1682c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
168318c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
16849ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
16857e6512fdSJed Brown   pc_gamg->use_sa_esteig    = PETSC_TRUE;
168618c3aa7eSMark   pc_gamg->emin             = 0;
168718c3aa7eSMark   pc_gamg->emax             = 0;
168818c3aa7eSMark 
1689c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
16909d5b6da9SMark F. Adams 
1691bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1692bd94a7aaSJed Brown   ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr);
16935b89ad90SMark F. Adams   PetscFunctionReturn(0);
16945b89ad90SMark F. Adams }
16953e3471ccSMark Adams 
16963e3471ccSMark Adams /*@C
16973e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
16988a690491SBarry Smith     from PCInitializePackage().
16993e3471ccSMark Adams 
17003e3471ccSMark Adams  Level: developer
17013e3471ccSMark Adams 
17023e3471ccSMark Adams  .seealso: PetscInitialize()
17033e3471ccSMark Adams @*/
17043e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
17053e3471ccSMark Adams {
17063e3471ccSMark Adams   PetscErrorCode ierr;
17074555aa8cSStefano Zampini   PetscInt       l;
17083e3471ccSMark Adams 
17093e3471ccSMark Adams   PetscFunctionBegin;
17103e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
17113e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
17123e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr);
17133e3471ccSMark Adams   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr);
17148e6d0c30SPeter Brune   ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr);
17153e3471ccSMark Adams   ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr);
1716c1c463dbSMark Adams 
1717c1c463dbSMark Adams   /* general events */
1718fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr);
1719fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr);
1720fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr);
1721fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr);
1722c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr);
1723c1c463dbSMark Adams   ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr);
1724fd1112cbSBarry Smith   ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr);
1725c1c463dbSMark Adams 
17265b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr);
172727436d3eSMark Adams   ierr = PetscLogEventRegister("  Create Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr);
172827436d3eSMark Adams   ierr = PetscLogEventRegister("  Filter Graph", PC_CLASSID, &petsc_gamg_setup_events[SET16]);CHKERRQ(ierr);
17295b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
17305b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
17315b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
17325b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr);
17335b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr);
17345b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr);
1735bb235841SBarry Smith   ierr = PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr);
17365b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr);
17375b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr);
17385b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr);
17395b89ad90SMark F. Adams   ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr);
17405b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr);
17415b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr);
17425b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr);
17435b89ad90SMark F. Adams   ierr = PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr);
17444555aa8cSStefano Zampini   for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
17454555aa8cSStefano Zampini     char ename[32];
17465b89ad90SMark F. Adams 
17474555aa8cSStefano Zampini     ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l);CHKERRQ(ierr);
17484555aa8cSStefano Zampini     ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);CHKERRQ(ierr);
17494555aa8cSStefano Zampini     ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l);CHKERRQ(ierr);
17504555aa8cSStefano Zampini     ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);CHKERRQ(ierr);
17514555aa8cSStefano Zampini     ierr = PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l);CHKERRQ(ierr);
17524555aa8cSStefano Zampini     ierr = PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);CHKERRQ(ierr);
17534555aa8cSStefano Zampini   }
17545b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
17555b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
17565b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
17575b89ad90SMark F. Adams   /* create timer stages */
17584555aa8cSStefano Zampini #if defined(GAMG_STAGES)
17595b89ad90SMark F. Adams   {
17605b89ad90SMark F. Adams     char     str[32];
17615b89ad90SMark F. Adams     PetscInt lidx;
17625b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
17635b89ad90SMark F. Adams     ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr);
17645b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
1765c9567895SMark       sprintf(str,"MG Level %d",(int)lidx);
17665b89ad90SMark F. Adams       ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr);
17675b89ad90SMark F. Adams     }
17685b89ad90SMark F. Adams   }
17695b89ad90SMark F. Adams #endif
17703e3471ccSMark Adams   PetscFunctionReturn(0);
17713e3471ccSMark Adams }
17723e3471ccSMark Adams 
17733e3471ccSMark Adams /*@C
17741c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
17751c1aac46SBarry Smith     called from PetscFinalize() automatically.
17763e3471ccSMark Adams 
17773e3471ccSMark Adams  Level: developer
17783e3471ccSMark Adams 
17793e3471ccSMark Adams  .seealso: PetscFinalize()
17803e3471ccSMark Adams @*/
17813e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
17823e3471ccSMark Adams {
17833e3471ccSMark Adams   PetscErrorCode ierr;
17843e3471ccSMark Adams 
17853e3471ccSMark Adams   PetscFunctionBegin;
17863e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
17873e3471ccSMark Adams   ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr);
17883e3471ccSMark Adams   PetscFunctionReturn(0);
17893e3471ccSMark Adams }
1790a36cf38bSToby Isaac 
1791a36cf38bSToby Isaac /*@C
1792a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1793a36cf38bSToby Isaac 
1794a36cf38bSToby Isaac  Input Parameters:
1795a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1796a36cf38bSToby Isaac  - create - function for creating the gamg context.
1797a36cf38bSToby Isaac 
1798a36cf38bSToby Isaac   Level: advanced
1799a36cf38bSToby Isaac 
18001c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1801a36cf38bSToby Isaac @*/
1802a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1803a36cf38bSToby Isaac {
1804a36cf38bSToby Isaac   PetscErrorCode ierr;
1805a36cf38bSToby Isaac 
1806a36cf38bSToby Isaac   PetscFunctionBegin;
1807a36cf38bSToby Isaac   ierr = PCGAMGInitializePackage();CHKERRQ(ierr);
1808a36cf38bSToby Isaac   ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr);
1809a36cf38bSToby Isaac   PetscFunctionReturn(0);
1810a36cf38bSToby Isaac }
1811a36cf38bSToby Isaac 
1812