xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 7827d75ba736e00c19d105c058b1c2ddcca945c7)
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 {
37d3d6bff4SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
38d3d6bff4SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
39d3d6bff4SMark F. Adams 
40d3d6bff4SMark F. Adams   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
421c1aac46SBarry Smith   pc_gamg->data_sz = 0;
439566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->orig_data));
445f80ce2aSJacob Faibussowitsch   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
4518c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
4618c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
4718c3aa7eSMark   }
4818c3aa7eSMark   pc_gamg->emin = 0;
4918c3aa7eSMark   pc_gamg->emax = 0;
50a2f3521dSMark F. Adams   PetscFunctionReturn(0);
51a2f3521dSMark F. Adams }
52a2f3521dSMark F. Adams 
535b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
545b89ad90SMark F. Adams /*
55c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
56a147abb0SMark F. Adams      of active processors.
575b89ad90SMark F. Adams 
585b89ad90SMark F. Adams    Input Parameter:
59a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
60a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
619d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
62c5bfad50SMark F. Adams    . cr_bs - coarse block size
633530afc2SMark F. Adams    In/Output Parameter:
64a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
65afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
6611e60469SMark F. Adams    Output Parameter:
673530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
685b89ad90SMark F. Adams */
695cb416c2SMark F. Adams 
70171cca9aSMark 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)
715b89ad90SMark F. Adams {
729d5b6da9SMark F. Adams   PC_MG           *mg         = (PC_MG*)pc->data;
73486a8d0bSJed Brown   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
74a2f3521dSMark F. Adams   Mat             Cmat,Pold=*a_P_inout;
753b4367a7SBarry Smith   MPI_Comm        comm;
76c5df96a5SBarry Smith   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
773ae0bb68SMark Adams   PetscInt        ncrs_eq,ncrs,f_bs;
785b89ad90SMark F. Adams 
795b89ad90SMark F. Adams   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine,&comm));
819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
839566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
859566063dSJacob Faibussowitsch   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
87038e3b61SMark F. Adams 
88ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
89ce7c7f2fSMark Adams 
903ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
919566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
923ae0bb68SMark Adams   if (pc_gamg->data_cell_rows>0) {
933ae0bb68SMark Adams     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
9473911c69SBarry Smith   } else {
953ae0bb68SMark Adams     PetscInt  bs;
969566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(Cmat, &bs));
973ae0bb68SMark Adams     ncrs = ncrs_eq/bs;
983ae0bb68SMark Adams   }
99c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
100c9567895SMark   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. */
101c9567895SMark #if defined(PETSC_HAVE_CUDA)
102c9567895SMark     PetscShmComm pshmcomm;
103c9567895SMark     PetscMPIInt  locrank;
104c9567895SMark     MPI_Comm     loccomm;
105c9567895SMark     PetscInt     s_nnodes,r_nnodes, new_new_size;
106c9567895SMark     cudaError_t  cerr;
107c9567895SMark     int          devCount;
1089566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGet(comm,&pshmcomm));
1099566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm,&loccomm));
1109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
111c9567895SMark     s_nnodes = !locrank;
1129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm));
113*7827d75bSBarry Smith     PetscCheck((size%r_nnodes) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes);
114c9567895SMark     devCount = 0;
115c9567895SMark     cerr = cudaGetDeviceCount(&devCount);
116c9567895SMark     cudaGetLastError(); /* Reset the last error */
117c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
118c9567895SMark       new_new_size = r_nnodes * devCount;
119c9567895SMark       new_size = new_new_size;
1209566063dSJacob Faibussowitsch       PetscCall(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));
121c9567895SMark     } else {
1229566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: With Cuda but no device. Use heuristics."));
123c9567895SMark       goto HEURISTIC;
124c9567895SMark     }
125c9567895SMark #else
126c9567895SMark     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
127c9567895SMark #endif
128c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
1292c71b3e2SJacob 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]);
130c9567895SMark     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
1319566063dSJacob Faibussowitsch     PetscCall(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]));
132c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
133c9567895SMark     new_size = 1;
1349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size));
135c9567895SMark   } else {
136472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
137c9567895SMark #if defined(PETSC_HAVE_CUDA)
138c9567895SMark     HEURISTIC:
139c9567895SMark #endif
1409566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
141a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
1427f66b68fSMark Adams     if (!new_size) new_size = 1; /* not likely, posible? */
143c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
1449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size));
145a2f3521dSMark F. Adams   }
1462e3501ffSMark Adams   if (new_size==nactive) {
147ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
148ce7c7f2fSMark Adams     if (new_size < size) {
149ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
1509566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive));
151ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
1529566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
1539566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
154ce7c7f2fSMark Adams       }
155ce7c7f2fSMark Adams     }
156ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1572e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
158192c0e8bSMark Adams     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
159885364a3SMark Adams     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
16071959b99SBarry Smith     nloc_old = ncrs_eq/cr_bs;
1612c71b3e2SJacob 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);
162ce7c7f2fSMark Adams     /* get new_size and rfactor */
163ce7c7f2fSMark Adams     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
164ce7c7f2fSMark Adams       /* find factor */
165ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
166ce7c7f2fSMark Adams       else {
167ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
168ce7c7f2fSMark Adams         jj = -1;
169ce7c7f2fSMark Adams         for (kk = 1 ; kk <= size ; kk++) {
170ce7c7f2fSMark Adams           if (!(size%kk)) { /* a candidate */
171ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
172ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
173ce7c7f2fSMark Adams             if (fact > best_fact) {
174ce7c7f2fSMark Adams               best_fact = fact; jj = kk;
175ce7c7f2fSMark Adams             }
176ce7c7f2fSMark Adams           }
177ce7c7f2fSMark Adams         }
178ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
179ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
180ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
181ce7c7f2fSMark Adams         else expand_factor = rfactor;
182ce7c7f2fSMark Adams       }
183ce7c7f2fSMark Adams       new_size = size/rfactor; /* make new size one that is factor */
1844cdfd227SMark       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1854cdfd227SMark         *a_Amat_crs = Cmat;
1869566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",((PetscObject)pc)->prefix,new_size,ncrs_eq));
187ce7c7f2fSMark Adams         PetscFunctionReturn(0);
188ce7c7f2fSMark Adams       }
189ce7c7f2fSMark Adams     }
1909566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0));
191a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
1929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &counts));
1932e3501ffSMark Adams     if (pc_gamg->repart) {
194a5b23f4aSJose E. Roman       /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
1955a9b9e01SMark F. Adams       Mat      adj;
1969566063dSJacob Faibussowitsch       PetscCall(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"));
197a2f3521dSMark F. Adams       /* get 'adj' */
198c5bfad50SMark F. Adams       if (cr_bs == 1) {
1999566063dSJacob Faibussowitsch         PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
200806fa848SBarry Smith       } else {
201a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
202eb07cef2SMark F. Adams         Mat               tMat;
203a2f3521dSMark F. Adams         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
204b4fbaa2aSMark F. Adams         const PetscScalar *vals;
205b4fbaa2aSMark F. Adams         const PetscInt    *idx;
206a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
20739d09545SMark Adams         static PetscInt   llev = 0; /* ugly but just used for debugging */
208d9558ea9SBarry Smith         MatType           mtype;
209b4fbaa2aSMark F. Adams 
2109566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz));
2119566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs));
2129566063dSJacob Faibussowitsch         PetscCall(MatGetSize(Cmat, &M, &N));
213c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
2149566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat,Ii,&ncols,NULL,NULL));
215c5bfad50SMark F. Adams           d_nnz[jj] = ncols/cr_bs;
216c5bfad50SMark F. Adams           o_nnz[jj] = ncols/cr_bs;
2179566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL));
2183ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
2193ae0bb68SMark Adams           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
22058471d46SMark F. Adams         }
2216876a03eSMark F. Adams 
2229566063dSJacob Faibussowitsch         PetscCall(MatGetType(Amat_fine,&mtype));
2239566063dSJacob Faibussowitsch         PetscCall(MatCreate(comm, &tMat));
2249566063dSJacob Faibussowitsch         PetscCall(MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE));
2259566063dSJacob Faibussowitsch         PetscCall(MatSetType(tMat,mtype));
2269566063dSJacob Faibussowitsch         PetscCall(MatSeqAIJSetPreallocation(tMat,0,d_nnz));
2279566063dSJacob Faibussowitsch         PetscCall(MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz));
2289566063dSJacob Faibussowitsch         PetscCall(PetscFree2(d_nnz,o_nnz));
229eb07cef2SMark F. Adams 
230a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
231c5bfad50SMark F. Adams           PetscInt dest_row = ii/cr_bs;
2329566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat,ii,&ncols,&idx,&vals));
233eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
234c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj]/cr_bs;
235eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
2369566063dSJacob Faibussowitsch             PetscCall(MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES));
237eb07cef2SMark F. Adams           }
2389566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat,ii,&ncols,&idx,&vals));
239eb07cef2SMark F. Adams         }
2409566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY));
2419566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY));
242eb07cef2SMark F. Adams 
243b4fbaa2aSMark F. Adams         if (llev++ == -1) {
244b4fbaa2aSMark F. Adams           PetscViewer viewer; char fname[32];
2459566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev));
2463b4367a7SBarry Smith           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
2479566063dSJacob Faibussowitsch           PetscCall(MatView(tMat, viewer));
2489566063dSJacob Faibussowitsch           PetscCall(PetscViewerDestroy(&viewer));
249b4fbaa2aSMark F. Adams         }
2509566063dSJacob Faibussowitsch         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
2519566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tMat));
252a2f3521dSMark F. Adams       } /* create 'adj' */
253f150b916SMark F. Adams 
254a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2555a9b9e01SMark F. Adams         char            prefix[256];
2565a9b9e01SMark F. Adams         const char      *pcpre;
257b4fbaa2aSMark F. Adams         const PetscInt  *is_idx;
258b4fbaa2aSMark F. Adams         MatPartitioning mpart;
259a4b7d37bSMark F. Adams         IS              proc_is;
2602f03bc48SMark F. Adams 
2619566063dSJacob Faibussowitsch         PetscCall(MatPartitioningCreate(comm, &mpart));
2629566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
2639566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
2649566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
2659566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix));
2669566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetFromOptions(mpart));
2679566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, new_size));
2689566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &proc_is));
2699566063dSJacob Faibussowitsch         PetscCall(MatPartitioningDestroy(&mpart));
2705a9b9e01SMark F. Adams 
2715ef31b24SMark F. Adams         /* collect IS info */
2729566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
2739566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(proc_is, &is_idx));
274a2f3521dSMark F. Adams         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
275c5bfad50SMark F. Adams           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
276ce7c7f2fSMark Adams             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
277eb07cef2SMark F. Adams           }
2785ef31b24SMark F. Adams         }
2799566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(proc_is, &is_idx));
2809566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&proc_is));
2815ef31b24SMark F. Adams       }
2829566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&adj));
2835a9b9e01SMark F. Adams 
2849566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
2859566063dSJacob Faibussowitsch       PetscCall(PetscFree(newproc_idx));
28631cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
287ce7c7f2fSMark Adams       PetscInt targetPE;
2882c71b3e2SJacob Faibussowitsch       PetscCheckFalse(new_size==nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
2899566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Number of equations (loc) %D with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq));
290ce7c7f2fSMark Adams       targetPE = (rank/rfactor)*expand_factor;
2919566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc));
292a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
293e33ef3b1SMark F. Adams 
29411e60469SMark F. Adams     /*
295a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
29611e60469SMark F. Adams     */
2979566063dSJacob Faibussowitsch     PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num));
2987700e67bSMark Adams     is_eq_num_prim = is_eq_num;
29911e60469SMark F. Adams     /*
300a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
30111e60469SMark F. Adams     */
3029566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(is_eq_newproc, size, counts));
303c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
3049566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_newproc));
305ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new/cr_bs;
306a2f3521dSMark F. Adams 
3079566063dSJacob Faibussowitsch     PetscCall(PetscFree(counts));
308885364a3SMark 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 */
309885364a3SMark Adams     {
310885364a3SMark Adams       Vec            src_crd, dest_crd;
311885364a3SMark 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;
312885364a3SMark Adams       VecScatter     vecscat;
313885364a3SMark Adams       PetscScalar    *array;
314885364a3SMark Adams       IS isscat;
315a2f3521dSMark F. Adams       /* move data (for primal equations only) */
31622063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
3179566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &dest_crd));
3189566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE));
3199566063dSJacob Faibussowitsch       PetscCall(VecSetType(dest_crd,VECSTANDARD)); /* this is needed! */
32011e60469SMark F. Adams       /*
3219d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
322c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
32311e60469SMark F. Adams       */
3249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(ncrs*node_data_sz, &tidx));
3259566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(is_eq_num_prim, &idx));
3263ae0bb68SMark Adams       for (ii=0,jj=0; ii<ncrs; ii++) {
327c5bfad50SMark F. Adams         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
328a2f3521dSMark F. Adams         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
32911e60469SMark F. Adams       }
3309566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(is_eq_num_prim, &idx));
3319566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat));
3329566063dSJacob Faibussowitsch       PetscCall(PetscFree(tidx));
33311e60469SMark F. Adams       /*
33411e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
33511e60469SMark F. Adams       */
3369566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd));
3379d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3383ae0bb68SMark Adams         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
3393ae0bb68SMark Adams         for (ii=0; ii<ncrs; ii++) {
3409d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
341a2f3521dSMark F. Adams             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
342c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
3439566063dSJacob Faibussowitsch             PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES));
344d3d6bff4SMark F. Adams           }
345038e3b61SMark F. Adams         }
346eb07cef2SMark F. Adams       }
3479566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(src_crd));
3489566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(src_crd));
34911e60469SMark F. Adams       /*
35011e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
35111e60469SMark F. Adams         to the correct processor
35211e60469SMark F. Adams       */
3539566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat));
3549566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isscat));
3559566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD));
3569566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD));
3579566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&vecscat));
3589566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&src_crd));
35911e60469SMark F. Adams       /*
36011e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
36111e60469SMark F. Adams       */
3629566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
3639566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data));
3642fa5cd67SKarl Rupp 
3653ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz*ncrs_new;
3663ae0bb68SMark Adams       strideNew        = ncrs_new*ndata_rows;
3672fa5cd67SKarl Rupp 
3689566063dSJacob Faibussowitsch       PetscCall(VecGetArray(dest_crd, &array));
3699d5b6da9SMark F. Adams       for (jj=0; jj<ndata_cols; jj++) {
3703ae0bb68SMark Adams         for (ii=0; ii<ncrs_new; ii++) {
3719d5b6da9SMark F. Adams           for (kk=0; kk<ndata_rows; kk++) {
372a2f3521dSMark F. Adams             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
373c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
374d3d6bff4SMark F. Adams           }
375038e3b61SMark F. Adams         }
376038e3b61SMark F. Adams       }
3779566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(dest_crd, &array));
3789566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&dest_crd));
379885364a3SMark Adams     }
380a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
3819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0));
38211e60469SMark F. Adams     /*
3837dae84e0SHong Zhang       Invert for MatCreateSubMatrix
38411e60469SMark F. Adams     */
3859566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
3869566063dSJacob Faibussowitsch     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
3879566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
388a2f3521dSMark F. Adams     if (is_eq_num != is_eq_num_prim) {
3899566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */
390a2f3521dSMark F. Adams     }
3913cb8563fSToby Isaac     if (Pcolumnperm) {
3929566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
3933cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3943cb8563fSToby Isaac     }
3959566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_num));
3969566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0));
3979566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0));
398a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
399a2f3521dSMark F. Adams     {
400a2f3521dSMark F. Adams       Mat       mat;
40190db8557SMark Adams       PetscBool flg;
4029566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
4039566063dSJacob Faibussowitsch       PetscCall(MatGetOption(Cmat, MAT_SPD, &flg));
40490db8557SMark Adams       if (flg) {
4059566063dSJacob Faibussowitsch         PetscCall(MatSetOption(mat, MAT_SPD,PETSC_TRUE));
40690db8557SMark Adams       } else {
4079566063dSJacob Faibussowitsch         PetscCall(MatGetOption(Cmat, MAT_HERMITIAN, &flg));
40890db8557SMark Adams         if (flg) {
4099566063dSJacob Faibussowitsch           PetscCall(MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE));
41090db8557SMark Adams         } else {
41190db8557SMark Adams #if !defined(PETSC_USE_COMPLEX)
4129566063dSJacob Faibussowitsch           PetscCall(MatGetOption(Cmat, MAT_SYMMETRIC, &flg));
41390db8557SMark Adams           if (flg) {
4149566063dSJacob Faibussowitsch             PetscCall(MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE));
41590db8557SMark Adams           }
41690db8557SMark Adams #endif
41790db8557SMark Adams         }
41890db8557SMark Adams       }
419a2f3521dSMark F. Adams       *a_Amat_crs = mat;
420a2f3521dSMark F. Adams     }
4219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cmat));
422a2f3521dSMark F. Adams 
4239566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0));
42411e60469SMark F. Adams     /* prolongator */
42511e60469SMark F. Adams     {
42611e60469SMark F. Adams       IS       findices;
427a2f3521dSMark F. Adams       PetscInt Istart,Iend;
428a2f3521dSMark F. Adams       Mat      Pnew;
42962294041SBarry Smith 
4309566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
4319566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0));
4329566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices));
4339566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(findices,f_bs));
4349566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
4359566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&findices));
4369566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
437c5bfad50SMark F. Adams 
4389566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0));
4399566063dSJacob Faibussowitsch       PetscCall(MatDestroy(a_P_inout));
440a2f3521dSMark F. Adams 
441a2f3521dSMark F. Adams       /* output - repartitioned */
442a2f3521dSMark F. Adams       *a_P_inout = Pnew;
443e33ef3b1SMark F. Adams     }
4449566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&new_eq_indices));
4455b89ad90SMark F. Adams 
446c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
447ce7c7f2fSMark Adams 
448ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
449ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
450ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4518bca76a6SMark Adams       static PetscInt llev = 2;
4529566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Pinning level %D to the CPU\n",((PetscObject)pc)->prefix,llev++));
453ce7c7f2fSMark Adams #endif
4549566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
4559566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
456adf5291fSStefano Zampini       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
457ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
458ce7c7f2fSMark Adams         PetscMPIInt size;
4599566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
460ce7c7f2fSMark Adams         if (size > 1) {
461ce7c7f2fSMark Adams           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
4629566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE));
4639566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE));
464ce7c7f2fSMark Adams         }
465ce7c7f2fSMark Adams       }
466ce7c7f2fSMark Adams     }
4679566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0));
468a2f3521dSMark F. Adams   }
4695b89ad90SMark F. Adams   PetscFunctionReturn(0);
4705b89ad90SMark F. Adams }
4715b89ad90SMark F. Adams 
4724b1575e2SStefano Zampini PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
4734b1575e2SStefano Zampini {
4744b1575e2SStefano Zampini   const char     *prefix;
4754b1575e2SStefano Zampini   char           addp[32];
4764b1575e2SStefano Zampini   PC_MG          *mg      = (PC_MG*)a_pc->data;
4774b1575e2SStefano Zampini   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4784b1575e2SStefano Zampini 
4794b1575e2SStefano Zampini   PetscFunctionBegin;
4809566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(a_pc,&prefix));
4819566063dSJacob Faibussowitsch   PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %D\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1));
4829566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2));
4839566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(*Gmat2,prefix));
4849566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level));
4859566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(*Gmat2,addp));
486b4da3a1bSStefano Zampini   if ((*Gmat2)->structurally_symmetric) {
4879566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB));
488b4da3a1bSStefano Zampini   } else {
4899566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
4909566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB));
491b4da3a1bSStefano Zampini   }
4929566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(*Gmat2));
4939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
4949566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(*Gmat2));
4959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
4969566063dSJacob Faibussowitsch   PetscCall(MatProductClear(*Gmat2));
4974b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
4984b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
4994b1575e2SStefano Zampini   PetscFunctionReturn(0);
5004b1575e2SStefano Zampini }
5014b1575e2SStefano Zampini 
5025b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
5035b89ad90SMark F. Adams /*
5045b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
5055b89ad90SMark F. Adams                     by setting data structures and options.
5065b89ad90SMark F. Adams 
5075b89ad90SMark F. Adams    Input Parameter:
5085b89ad90SMark F. Adams .  pc - the preconditioner context
5095b89ad90SMark F. Adams 
5105b89ad90SMark F. Adams */
5119d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc)
5125b89ad90SMark F. Adams {
5135b89ad90SMark F. Adams   PetscErrorCode ierr;
5149d5b6da9SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
5155b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
5162adcac29SMark F. Adams   Mat            Pmat     = pc->pmat;
51718c3aa7eSMark   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
5183b4367a7SBarry Smith   MPI_Comm       comm;
519c5df96a5SBarry Smith   PetscMPIInt    rank,size,nactivepe;
52018c3aa7eSMark   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
52118c3aa7eSMark   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
522a2f3521dSMark F. Adams   PetscLogDouble nnz0=0.,nnztot=0.;
523569f4572SMark Adams   MatInfo        info;
524171cca9aSMark Adams   PetscBool      is_last = PETSC_FALSE;
5255ef31b24SMark F. Adams 
5265b89ad90SMark F. Adams   PetscFunctionBegin;
5279566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
5289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
5299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
530dfd5c07aSMark F. Adams 
5318abdc6daSStefano Zampini   if (pc->setupcalled) {
5328abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
533878e152fSMark F. Adams       /* reset everything */
5349566063dSJacob Faibussowitsch       PetscCall(PCReset_MG(pc));
535878e152fSMark F. Adams       pc->setupcalled = 0;
536806fa848SBarry Smith     } else {
53784d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
53803a628feSMark F. Adams       /* just do Galerkin grids */
53958471d46SMark F. Adams       Mat          B,dA,dB;
54058471d46SMark F. Adams 
5419d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
5424555aa8cSStefano Zampini         PetscInt gl;
54358471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5449566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB));
54558471d46SMark F. Adams         /* (re)set to get dirty flag */
5469566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB));
54758471d46SMark F. Adams 
5484555aa8cSStefano Zampini         for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
5498abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX ;
5508abdc6daSStefano Zampini 
5518abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
5529566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B));
5538abdc6daSStefano Zampini           if (B->product) {
5548abdc6daSStefano Zampini             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
5558abdc6daSStefano Zampini               reuse = MAT_REUSE_MATRIX;
55603a628feSMark F. Adams             }
5578abdc6daSStefano Zampini           }
5589566063dSJacob Faibussowitsch           if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A));
5598abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
5609566063dSJacob Faibussowitsch             PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %D\n",((PetscObject)pc)->prefix,level));
5618abdc6daSStefano Zampini           } else {
5629566063dSJacob Faibussowitsch             PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %D\n",((PetscObject)pc)->prefix,level));
5638abdc6daSStefano Zampini           }
5649566063dSJacob Faibussowitsch           PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
5659566063dSJacob Faibussowitsch           PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B));
5669566063dSJacob Faibussowitsch           PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
56763b77682SMark Adams           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
5689566063dSJacob Faibussowitsch           PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B));
56958471d46SMark F. Adams           dB   = B;
57058471d46SMark F. Adams         }
5715f8cf99dSMark F. Adams       }
572d5280255SMark F. Adams 
5739566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
57458471d46SMark F. Adams       PetscFunctionReturn(0);
575eb07cef2SMark F. Adams     }
576878e152fSMark F. Adams   }
577f6536408SMark F. Adams 
578878e152fSMark F. Adams   if (!pc_gamg->data) {
579878e152fSMark F. Adams     if (pc_gamg->orig_data) {
5809566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(Pmat, &bs));
5819566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
5822fa5cd67SKarl Rupp 
583878e152fSMark F. Adams       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
584878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
585878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5862fa5cd67SKarl Rupp 
5879566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
588878e152fSMark F. Adams       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
589806fa848SBarry Smith     } else {
5905f80ce2aSJacob Faibussowitsch       PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
5919566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat));
5929d5b6da9SMark F. Adams     }
593878e152fSMark F. Adams   }
594878e152fSMark F. Adams 
595878e152fSMark F. Adams   /* cache original data for reuse */
5961c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
5979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
598878e152fSMark F. Adams     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
599878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
600878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
601878e152fSMark F. Adams   }
602038e3b61SMark F. Adams 
603302f38e8SMark F. Adams   /* get basic dims */
6049566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Pmat, &bs));
6059566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Pmat, &M, &N));
60684d3f75bSMark F. Adams 
6079566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */
608569f4572SMark Adams   nnz0   = info.nz_used;
609569f4572SMark Adams   nnztot = info.nz_used;
6109566063dSJacob Faibussowitsch   PetscCall(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));
611569f4572SMark Adams 
612a2f3521dSMark F. Adams   /* Get A_i and R_i */
61362294041SBarry Smith   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
6149ab59c8bSMark Adams     pc_gamg->current_level = level;
6155f80ce2aSJacob Faibussowitsch     PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
6165b89ad90SMark F. Adams     level1 = level + 1;
6179566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0));
6184555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6199566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(gamg_stages[level]));
620b4fbaa2aSMark F. Adams #endif
621c8b0795cSMark F. Adams     { /* construct prolongator */
622725b86d8SJed Brown       Mat              Gmat;
6230cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6247700e67bSMark Adams       Mat              Prol11;
625c8b0795cSMark F. Adams 
6269566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat));
6279566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists));
6289566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11));
629c8b0795cSMark F. Adams 
630a2f3521dSMark F. Adams       /* could have failed to create new level */
631a2f3521dSMark F. Adams       if (Prol11) {
632f7df55f0SStefano Zampini         const char *prefix;
633f7df55f0SStefano Zampini         char       addp[32];
634f7df55f0SStefano Zampini 
6359d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6369566063dSJacob Faibussowitsch         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs));
637a2f3521dSMark F. Adams 
638fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
639c8b0795cSMark F. Adams           /* smooth */
6409566063dSJacob Faibussowitsch           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
641c8b0795cSMark F. Adams         }
642c8b0795cSMark F. Adams 
6430c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
6441b18a24aSMark Adams           PetscInt bs;
6459566063dSJacob Faibussowitsch           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL));
6469566063dSJacob Faibussowitsch           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
647ffc955d6SMark F. Adams         }
648ffc955d6SMark F. Adams 
6499566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc,&prefix));
6509566063dSJacob Faibussowitsch         PetscCall(MatSetOptionsPrefix(Prol11,prefix));
6519566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level));
6529566063dSJacob Faibussowitsch         PetscCall(MatAppendOptionsPrefix(Prol11,addp));
65391f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
654f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
6559566063dSJacob Faibussowitsch         PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
6569566063dSJacob Faibussowitsch         PetscCall(MatSetFromOptions(Prol11));
6574bde40a0SMark Adams         Parr[level1] = Prol11;
6584bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6594bde40a0SMark Adams 
6609566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat));
6619566063dSJacob Faibussowitsch       PetscCall(PetscCDDestroy(agg_lists));
662a2f3521dSMark F. Adams     } /* construct prolongator scope */
6639566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0));
6647f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
665171cca9aSMark Adams     if (!Parr[level1]) { /* failed to coarsen */
6669566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Stop gridding, level %D\n",((PetscObject)pc)->prefix,level));
6674555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6689566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePop());
669a90e85d9SMark Adams #endif
670c8b0795cSMark F. Adams       break;
671c8b0795cSMark F. Adams     }
6729566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0));
6739566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
6745f80ce2aSJacob Faibussowitsch     PetscCheck(!is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
675171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6760e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
6779566063dSJacob Faibussowitsch     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
678a2f3521dSMark F. Adams 
6799566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0));
6809566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
6819566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
682569f4572SMark Adams     nnztot += info.nz_used;
6839566063dSJacob Faibussowitsch     PetscCall(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));
684569f4572SMark Adams 
6854555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6869566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
687b4fbaa2aSMark F. Adams #endif
688a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
6899ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
6909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",((PetscObject)pc)->prefix,level,M/bs));
691a90e85d9SMark Adams       level++;
692a90e85d9SMark Adams       break;
693a90e85d9SMark Adams     }
694c8b0795cSMark F. Adams   } /* levels */
6959566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
696c8b0795cSMark F. Adams 
6979566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"%s: %D levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0));
6989d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
6995b89ad90SMark F. Adams   fine_level       = level;
7009566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL));
7015b89ad90SMark F. Adams 
70262294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7030ed2132dSStefano Zampini 
704d5280255SMark F. Adams     /* set default smoothers & set operators */
70562294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
706ffc955d6SMark F. Adams       KSP smoother;
707ffc955d6SMark F. Adams       PC  subpc;
708a2f3521dSMark F. Adams 
7099566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7109566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(smoother, &subpc));
711ffc955d6SMark F. Adams 
7129566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
713a2f3521dSMark F. Adams       /* set ops */
7149566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
7159566063dSJacob Faibussowitsch       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1]));
716a2f3521dSMark F. Adams 
717a2f3521dSMark F. Adams       /* set defaults */
7189566063dSJacob Faibussowitsch       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
719a2f3521dSMark F. Adams 
7200c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
7210c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
7222d3561bbSSatish Balay         PetscInt sz;
7237a28f3e5SMark Adams         IS       *iss;
724a2f3521dSMark F. Adams 
7252d3561bbSSatish Balay         sz   = nASMBlocksArr[level];
7267a28f3e5SMark Adams         iss  = ASMLocalIDsArr[level];
7279566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCASM));
7289566063dSJacob Faibussowitsch         PetscCall(PCASMSetOverlap(subpc, 0));
7299566063dSJacob Faibussowitsch         PetscCall(PCASMSetType(subpc,PC_ASM_BASIC));
7307f66b68fSMark Adams         if (!sz) {
731ffc955d6SMark F. Adams           IS       is;
7329566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
7339566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
7349566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
735806fa848SBarry Smith         } else {
736a94c3b12SMark F. Adams           PetscInt kk;
7379566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss));
738a94c3b12SMark F. Adams           for (kk=0; kk<sz; kk++) {
7399566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&iss[kk]));
740a94c3b12SMark F. Adams           }
7419566063dSJacob Faibussowitsch           PetscCall(PetscFree(iss));
742ffc955d6SMark F. Adams         }
7430298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
744ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
745806fa848SBarry Smith       } else {
7469566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCJACOBI));
747ffc955d6SMark F. Adams       }
748d5280255SMark F. Adams     }
749d5280255SMark F. Adams     {
750d5280255SMark F. Adams       /* coarse grid */
751d5280255SMark F. Adams       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
752d5280255SMark F. Adams       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
7530ed2132dSStefano Zampini 
7549566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7559566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
756cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
7579566063dSJacob Faibussowitsch         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
7589566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(smoother, &subpc));
7599566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCBJACOBI));
7609566063dSJacob Faibussowitsch         PetscCall(PCSetUp(subpc));
7619566063dSJacob Faibussowitsch         PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2));
7625f80ce2aSJacob Faibussowitsch         PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
7639566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(k2[0],&pc2));
7649566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc2, PCLU));
7659566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS));
7669566063dSJacob Faibussowitsch         PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1));
7679566063dSJacob Faibussowitsch         PetscCall(KSPSetType(k2[0], KSPPREONLY));
768d5280255SMark F. Adams       }
769cf8ae1d3SMark Adams     }
770d5280255SMark F. Adams 
771d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
7729566063dSJacob Faibussowitsch     ierr = PetscObjectOptionsBegin((PetscObject)pc);PetscCall(ierr);
7739566063dSJacob Faibussowitsch     PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc));
7749566063dSJacob Faibussowitsch     ierr = PetscOptionsEnd();PetscCall(ierr);
7759566063dSJacob Faibussowitsch     PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
776d5280255SMark F. Adams 
77718c3aa7eSMark     /* setup cheby eigen estimates from SA */
7787e6512fdSJed Brown     if (pc_gamg->use_sa_esteig) {
77918c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
78018c3aa7eSMark         KSP       smoother;
78118c3aa7eSMark         PetscBool ischeb;
7820ed2132dSStefano Zampini 
7839566063dSJacob Faibussowitsch         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7849566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb));
78518c3aa7eSMark         if (ischeb) {
78618c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
7870ed2132dSStefano Zampini 
7882de708cbSJed Brown           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
789efe053fcSJed Brown           if (mg->max_eigen_DinvA[level] > 0) {
7907e6512fdSJed Brown             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
7917e6512fdSJed Brown             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
79218c3aa7eSMark             PetscReal emax,emin;
7930ed2132dSStefano Zampini 
79418c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
79518c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
7969566063dSJacob Faibussowitsch             PetscCall(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));
7970a94a983SJed Brown             cheb->emin_provided = emin;
7980a94a983SJed Brown             cheb->emax_provided = emax;
79918c3aa7eSMark           }
80018c3aa7eSMark         }
80118c3aa7eSMark       }
80218c3aa7eSMark     }
8030ed2132dSStefano Zampini 
8049566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8050ed2132dSStefano Zampini 
806d5280255SMark F. Adams     /* clean up */
807d5280255SMark F. Adams     for (level=1; level<pc_gamg->Nlevels; level++) {
8089566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Parr[level]));
8099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Aarr[level]));
8105b89ad90SMark F. Adams     }
811806fa848SBarry Smith   } else {
8125f8cf99dSMark F. Adams     KSP smoother;
8130ed2132dSStefano Zampini 
8149566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n"));
8159566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
8169566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
8179566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPPREONLY));
8189566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8195f8cf99dSMark F. Adams   }
8205b89ad90SMark F. Adams   PetscFunctionReturn(0);
8215b89ad90SMark F. Adams }
8225b89ad90SMark F. Adams 
823eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */
8245b89ad90SMark F. Adams /*
8255b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8265b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8275b89ad90SMark F. Adams 
8285b89ad90SMark F. Adams    Input Parameter:
8295b89ad90SMark F. Adams .  pc - the preconditioner context
8305b89ad90SMark F. Adams 
8315b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8325b89ad90SMark F. Adams */
8335b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
8345b89ad90SMark F. Adams {
8355b89ad90SMark F. Adams   PC_MG          *mg     = (PC_MG*)pc->data;
8365b89ad90SMark F. Adams   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
8375b89ad90SMark F. Adams 
8385b89ad90SMark F. Adams   PetscFunctionBegin;
8399566063dSJacob Faibussowitsch   PetscCall(PCReset_GAMG(pc));
8409b8ffb57SJed Brown   if (pc_gamg->ops->destroy) {
8419566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
8429b8ffb57SJed Brown   }
8439566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->ops));
8449566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
8459566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg));
8469566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
8475b89ad90SMark F. Adams   PetscFunctionReturn(0);
8485b89ad90SMark F. Adams }
8495b89ad90SMark F. Adams 
850676e1743SMark F. Adams /*@
851cab9ed1eSBarry Smith    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
852676e1743SMark F. Adams 
8531cc46a46SBarry Smith    Logically Collective on PC
854676e1743SMark F. Adams 
855676e1743SMark F. Adams    Input Parameters:
8561cc46a46SBarry Smith +  pc - the preconditioner context
8571cc46a46SBarry Smith -  n - the number of equations
858676e1743SMark F. Adams 
859676e1743SMark F. Adams    Options Database Key:
860147403d9SBarry Smith .  -pc_gamg_process_eq_limit <limit> - set the limit
861676e1743SMark F. Adams 
86295452b02SPatrick Sanan    Notes:
86395452b02SPatrick 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
864cab9ed1eSBarry Smith     that has degrees of freedom
865cab9ed1eSBarry Smith 
866676e1743SMark F. Adams    Level: intermediate
867676e1743SMark F. Adams 
868c9567895SMark .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
869676e1743SMark F. Adams @*/
870676e1743SMark F. Adams PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
871676e1743SMark F. Adams {
872676e1743SMark F. Adams   PetscFunctionBegin;
873676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8749566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n)));
875676e1743SMark F. Adams   PetscFunctionReturn(0);
876676e1743SMark F. Adams }
877676e1743SMark F. Adams 
8781e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
879676e1743SMark F. Adams {
880c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
881c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
882676e1743SMark F. Adams 
883676e1743SMark F. Adams   PetscFunctionBegin;
8849d5b6da9SMark F. Adams   if (n>0) pc_gamg->min_eq_proc = n;
885676e1743SMark F. Adams   PetscFunctionReturn(0);
886676e1743SMark F. Adams }
887676e1743SMark F. Adams 
888389730f3SMark F. Adams /*@
889cab9ed1eSBarry Smith    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
890389730f3SMark F. Adams 
891389730f3SMark F. Adams  Collective on PC
892389730f3SMark F. Adams 
893389730f3SMark F. Adams    Input Parameters:
8941cc46a46SBarry Smith +  pc - the preconditioner context
8951cc46a46SBarry Smith -  n - maximum number of equations to aim for
896389730f3SMark F. Adams 
897389730f3SMark F. Adams    Options Database Key:
898147403d9SBarry Smith .  -pc_gamg_coarse_eq_limit <limit> - set the limit
899389730f3SMark F. Adams 
90074329af1SBarry Smith    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
90174329af1SBarry Smith      has less than 1000 unknowns.
90274329af1SBarry Smith 
903389730f3SMark F. Adams    Level: intermediate
904389730f3SMark F. Adams 
905c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
906389730f3SMark F. Adams @*/
907389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
908389730f3SMark F. Adams {
909389730f3SMark F. Adams   PetscFunctionBegin;
910389730f3SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9119566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n)));
912389730f3SMark F. Adams   PetscFunctionReturn(0);
913389730f3SMark F. Adams }
914389730f3SMark F. Adams 
9151e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
916389730f3SMark F. Adams {
917389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
918389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
919389730f3SMark F. Adams 
920389730f3SMark F. Adams   PetscFunctionBegin;
9219d5b6da9SMark F. Adams   if (n>0) pc_gamg->coarse_eq_limit = n;
922389730f3SMark F. Adams   PetscFunctionReturn(0);
923389730f3SMark F. Adams }
924389730f3SMark F. Adams 
925676e1743SMark F. Adams /*@
926cab9ed1eSBarry Smith    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
927676e1743SMark F. Adams 
928676e1743SMark F. Adams    Collective on PC
929676e1743SMark F. Adams 
930676e1743SMark F. Adams    Input Parameters:
9311cc46a46SBarry Smith +  pc - the preconditioner context
9321cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
933676e1743SMark F. Adams 
934676e1743SMark F. Adams    Options Database Key:
935147403d9SBarry Smith .  -pc_gamg_repartition <true,false> - turn on the repartitioning
936676e1743SMark F. Adams 
93795452b02SPatrick Sanan    Notes:
93895452b02SPatrick Sanan     this will generally improve the loading balancing of the work on each level
939cab9ed1eSBarry Smith 
940676e1743SMark F. Adams    Level: intermediate
941676e1743SMark F. Adams 
942676e1743SMark F. Adams @*/
943cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
944676e1743SMark F. Adams {
945676e1743SMark F. Adams   PetscFunctionBegin;
946676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9479566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n)));
948676e1743SMark F. Adams   PetscFunctionReturn(0);
949676e1743SMark F. Adams }
950676e1743SMark F. Adams 
951cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
952676e1743SMark F. Adams {
953c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
954c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
955676e1743SMark F. Adams 
956676e1743SMark F. Adams   PetscFunctionBegin;
9579d5b6da9SMark F. Adams   pc_gamg->repart = n;
958676e1743SMark F. Adams   PetscFunctionReturn(0);
959676e1743SMark F. Adams }
960676e1743SMark F. Adams 
961dfd5c07aSMark F. Adams /*@
9627e6512fdSJed Brown    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother
96318c3aa7eSMark 
96418c3aa7eSMark    Collective on PC
96518c3aa7eSMark 
96618c3aa7eSMark    Input Parameters:
96718c3aa7eSMark +  pc - the preconditioner context
96818c3aa7eSMark -  n - number of its
96918c3aa7eSMark 
97018c3aa7eSMark    Options Database Key:
971147403d9SBarry Smith .  -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
97218c3aa7eSMark 
97318c3aa7eSMark    Notes:
9747e6512fdSJed 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$.
9757e6512fdSJed Brown    Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
9767e6512fdSJed Brown    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused.
977efe053fcSJed Brown    This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used.
978efe053fcSJed Brown    It became default in PETSc 3.17.
97918c3aa7eSMark 
9807e6512fdSJed Brown    Level: advanced
98118c3aa7eSMark 
98273f7197eSJed Brown .seealso: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet()
98318c3aa7eSMark @*/
98418c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
98518c3aa7eSMark {
98618c3aa7eSMark   PetscFunctionBegin;
98718c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9889566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n)));
98918c3aa7eSMark   PetscFunctionReturn(0);
99018c3aa7eSMark }
99118c3aa7eSMark 
9920ed2132dSStefano Zampini static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
99318c3aa7eSMark {
99418c3aa7eSMark   PC_MG   *mg      = (PC_MG*)pc->data;
99518c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
99618c3aa7eSMark 
99718c3aa7eSMark   PetscFunctionBegin;
9987e6512fdSJed Brown   pc_gamg->use_sa_esteig = n;
99918c3aa7eSMark   PetscFunctionReturn(0);
100018c3aa7eSMark }
100118c3aa7eSMark 
100218c3aa7eSMark /*@
100318c3aa7eSMark    PCGAMGSetEigenvalues - Set eigenvalues
100418c3aa7eSMark 
100518c3aa7eSMark    Collective on PC
100618c3aa7eSMark 
100718c3aa7eSMark    Input Parameters:
100818c3aa7eSMark +  pc - the preconditioner context
100918c3aa7eSMark -  emax - max eigenvalue
101018c3aa7eSMark .  emin - min eigenvalue
101118c3aa7eSMark 
101218c3aa7eSMark    Options Database Key:
1013147403d9SBarry Smith .  -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
101418c3aa7eSMark 
101518c3aa7eSMark    Level: intermediate
101618c3aa7eSMark 
101773f7197eSJed Brown .seealso: PCGAMGSetUseSAEstEig()
101818c3aa7eSMark @*/
101918c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
102018c3aa7eSMark {
102118c3aa7eSMark   PetscFunctionBegin;
102218c3aa7eSMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10239566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin)));
102418c3aa7eSMark   PetscFunctionReturn(0);
102518c3aa7eSMark }
102641ffd417SStefano Zampini 
102718c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
102818c3aa7eSMark {
102918c3aa7eSMark   PC_MG          *mg      = (PC_MG*)pc->data;
103018c3aa7eSMark   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
103118c3aa7eSMark 
103218c3aa7eSMark   PetscFunctionBegin;
10335f80ce2aSJacob Faibussowitsch   PetscCheck(emax > emin,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
10345f80ce2aSJacob Faibussowitsch   PetscCheck(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);
103518c3aa7eSMark   pc_gamg->emax = emax;
103618c3aa7eSMark   pc_gamg->emin = emin;
103718c3aa7eSMark   PetscFunctionReturn(0);
103818c3aa7eSMark }
103918c3aa7eSMark 
104018c3aa7eSMark /*@
1041cab9ed1eSBarry Smith    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1042dfd5c07aSMark F. Adams 
1043dfd5c07aSMark F. Adams    Collective on PC
1044dfd5c07aSMark F. Adams 
1045dfd5c07aSMark F. Adams    Input Parameters:
10461cc46a46SBarry Smith +  pc - the preconditioner context
10471cc46a46SBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
1048dfd5c07aSMark F. Adams 
1049dfd5c07aSMark F. Adams    Options Database Key:
1050147403d9SBarry Smith .  -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1051dfd5c07aSMark F. Adams 
1052dfd5c07aSMark F. Adams    Level: intermediate
1053dfd5c07aSMark F. Adams 
105495452b02SPatrick Sanan    Notes:
1055147403d9SBarry Smith     May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1056cab9ed1eSBarry Smith     rebuilding the preconditioner quicker.
1057cab9ed1eSBarry Smith 
1058dfd5c07aSMark F. Adams @*/
10591cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1060dfd5c07aSMark F. Adams {
1061dfd5c07aSMark F. Adams   PetscFunctionBegin;
1062dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10639566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n)));
1064dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1065dfd5c07aSMark F. Adams }
1066dfd5c07aSMark F. Adams 
10671cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1068dfd5c07aSMark F. Adams {
1069dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1070dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1071dfd5c07aSMark F. Adams 
1072dfd5c07aSMark F. Adams   PetscFunctionBegin;
1073dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
1074dfd5c07aSMark F. Adams   PetscFunctionReturn(0);
1075dfd5c07aSMark F. Adams }
1076dfd5c07aSMark F. Adams 
1077ffc955d6SMark F. Adams /*@
1078cab9ed1eSBarry 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.
1079ffc955d6SMark F. Adams 
1080ffc955d6SMark F. Adams    Collective on PC
1081ffc955d6SMark F. Adams 
1082ffc955d6SMark F. Adams    Input Parameters:
1083cab9ed1eSBarry Smith +  pc - the preconditioner context
1084cab9ed1eSBarry Smith -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1085ffc955d6SMark F. Adams 
1086ffc955d6SMark F. Adams    Options Database Key:
1087147403d9SBarry Smith .  -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1088ffc955d6SMark F. Adams 
1089ffc955d6SMark F. Adams    Level: intermediate
1090ffc955d6SMark F. Adams 
1091ffc955d6SMark F. Adams @*/
1092cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1093ffc955d6SMark F. Adams {
1094ffc955d6SMark F. Adams   PetscFunctionBegin;
1095ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10969566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg)));
1097ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1098ffc955d6SMark F. Adams }
1099ffc955d6SMark F. Adams 
1100cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1101ffc955d6SMark F. Adams {
1102ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1103ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1104ffc955d6SMark F. Adams 
1105ffc955d6SMark F. Adams   PetscFunctionBegin;
1106cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
1107ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1108ffc955d6SMark F. Adams }
1109ffc955d6SMark F. Adams 
1110171cca9aSMark Adams /*@
1111cf8ae1d3SMark Adams    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1112171cca9aSMark Adams 
1113171cca9aSMark Adams    Collective on PC
1114171cca9aSMark Adams 
1115171cca9aSMark Adams    Input Parameters:
1116171cca9aSMark Adams +  pc - the preconditioner context
1117cf8ae1d3SMark Adams -  flg - PETSC_TRUE to not force coarse grid onto one processor
1118171cca9aSMark Adams 
1119171cca9aSMark Adams    Options Database Key:
1120147403d9SBarry Smith .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1121171cca9aSMark Adams 
1122171cca9aSMark Adams    Level: intermediate
1123171cca9aSMark Adams 
112439d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1125171cca9aSMark Adams @*/
1126171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1127171cca9aSMark Adams {
1128171cca9aSMark Adams   PetscFunctionBegin;
1129171cca9aSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11309566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg)));
1131171cca9aSMark Adams   PetscFunctionReturn(0);
1132171cca9aSMark Adams }
1133171cca9aSMark Adams 
1134171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1135171cca9aSMark Adams {
1136171cca9aSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1137171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1138171cca9aSMark Adams 
1139171cca9aSMark Adams   PetscFunctionBegin;
1140171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
1141ffc955d6SMark F. Adams   PetscFunctionReturn(0);
1142ffc955d6SMark F. Adams }
1143ffc955d6SMark F. Adams 
11444ef23d27SMark F. Adams /*@
1145ce7c7f2fSMark Adams    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1146ce7c7f2fSMark Adams 
1147ce7c7f2fSMark Adams    Collective on PC
1148ce7c7f2fSMark Adams 
1149ce7c7f2fSMark Adams    Input Parameters:
1150ce7c7f2fSMark Adams +  pc - the preconditioner context
1151ce7c7f2fSMark Adams -  flg - PETSC_TRUE to pin coarse grids to CPU
1152ce7c7f2fSMark Adams 
1153ce7c7f2fSMark Adams    Options Database Key:
1154147403d9SBarry Smith .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1155ce7c7f2fSMark Adams 
1156ce7c7f2fSMark Adams    Level: intermediate
1157ce7c7f2fSMark Adams 
115839d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1159ce7c7f2fSMark Adams @*/
1160ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1161ce7c7f2fSMark Adams {
1162ce7c7f2fSMark Adams   PetscFunctionBegin;
1163ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11649566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg)));
1165ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1166ce7c7f2fSMark Adams }
1167ce7c7f2fSMark Adams 
1168ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1169ce7c7f2fSMark Adams {
1170ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1171ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1172ce7c7f2fSMark Adams 
1173ce7c7f2fSMark Adams   PetscFunctionBegin;
1174ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
1175ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1176ce7c7f2fSMark Adams }
1177ce7c7f2fSMark Adams 
1178ce7c7f2fSMark Adams /*@
1179147403d9SBarry Smith    PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1180ce7c7f2fSMark Adams 
1181ce7c7f2fSMark Adams    Collective on PC
1182ce7c7f2fSMark Adams 
1183ce7c7f2fSMark Adams    Input Parameters:
1184ce7c7f2fSMark Adams +  pc - the preconditioner context
1185ce7c7f2fSMark Adams -  flg - Layout type
1186ce7c7f2fSMark Adams 
1187ce7c7f2fSMark Adams    Options Database Key:
1188147403d9SBarry Smith .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1189ce7c7f2fSMark Adams 
1190ce7c7f2fSMark Adams    Level: intermediate
1191ce7c7f2fSMark Adams 
119239d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1193ce7c7f2fSMark Adams @*/
1194ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1195ce7c7f2fSMark Adams {
1196ce7c7f2fSMark Adams   PetscFunctionBegin;
1197ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11989566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg)));
1199ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1200ce7c7f2fSMark Adams }
1201ce7c7f2fSMark Adams 
1202ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1203ce7c7f2fSMark Adams {
1204ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1205ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1206ce7c7f2fSMark Adams 
1207ce7c7f2fSMark Adams   PetscFunctionBegin;
1208ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
1209ce7c7f2fSMark Adams   PetscFunctionReturn(0);
1210ce7c7f2fSMark Adams }
1211ce7c7f2fSMark Adams 
1212ce7c7f2fSMark Adams /*@
12131cc46a46SBarry Smith    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
12144ef23d27SMark F. Adams 
12154ef23d27SMark F. Adams    Not collective on PC
12164ef23d27SMark F. Adams 
12174ef23d27SMark F. Adams    Input Parameters:
12181cc46a46SBarry Smith +  pc - the preconditioner
12191cc46a46SBarry Smith -  n - the maximum number of levels to use
12204ef23d27SMark F. Adams 
12214ef23d27SMark F. Adams    Options Database Key:
1222147403d9SBarry Smith .  -pc_mg_levels <n> - set the maximum number of levels to allow
12234ef23d27SMark F. Adams 
12244ef23d27SMark F. Adams    Level: intermediate
12254ef23d27SMark F. Adams 
12264ef23d27SMark F. Adams @*/
12274ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
12284ef23d27SMark F. Adams {
12294ef23d27SMark F. Adams   PetscFunctionBegin;
12304ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12319566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n)));
12324ef23d27SMark F. Adams   PetscFunctionReturn(0);
12334ef23d27SMark F. Adams }
12344ef23d27SMark F. Adams 
12351e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
12364ef23d27SMark F. Adams {
12374ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
12384ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
12394ef23d27SMark F. Adams 
12404ef23d27SMark F. Adams   PetscFunctionBegin;
12419d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
12424ef23d27SMark F. Adams   PetscFunctionReturn(0);
12434ef23d27SMark F. Adams }
12444ef23d27SMark F. Adams 
12453542efc5SMark F. Adams /*@
12463542efc5SMark F. Adams    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
12473542efc5SMark F. Adams 
12483542efc5SMark F. Adams    Not collective on PC
12493542efc5SMark F. Adams 
12503542efc5SMark F. Adams    Input Parameters:
12511cc46a46SBarry Smith +  pc - the preconditioner context
1252c9567895SMark .  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
1253055c8bd0SJed Brown -  n - number of threshold values provided in array
12543542efc5SMark F. Adams 
12553542efc5SMark F. Adams    Options Database Key:
1256147403d9SBarry Smith .  -pc_gamg_threshold <threshold> - the threshold to drop edges
12573542efc5SMark F. Adams 
125895452b02SPatrick Sanan    Notes:
1259af3c827dSMark 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.
1260af3c827dSMark 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.
1261cab9ed1eSBarry Smith 
1262055c8bd0SJed Brown     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1263055c8bd0SJed Brown     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1264055c8bd0SJed Brown     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1265055c8bd0SJed Brown 
12663542efc5SMark F. Adams    Level: intermediate
12673542efc5SMark F. Adams 
1268af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
12693542efc5SMark F. Adams @*/
1270c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
12713542efc5SMark F. Adams {
12723542efc5SMark F. Adams   PetscFunctionBegin;
12733542efc5SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1274055c8bd0SJed Brown   if (n) PetscValidRealPointer(v,2);
12759566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n)));
12763542efc5SMark F. Adams   PetscFunctionReturn(0);
12773542efc5SMark F. Adams }
12783542efc5SMark F. Adams 
1279c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
12803542efc5SMark F. Adams {
1281c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1282c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1283c1eae691SMark Adams   PetscInt i;
1284c1eae691SMark Adams   PetscFunctionBegin;
1285055c8bd0SJed Brown   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1286055c8bd0SJed Brown   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1287c1eae691SMark Adams   PetscFunctionReturn(0);
1288c1eae691SMark Adams }
1289c1eae691SMark Adams 
1290c1eae691SMark Adams /*@
1291147403d9SBarry Smith    PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids
1292c9567895SMark 
1293c9567895SMark    Collective on PC
1294c9567895SMark 
1295c9567895SMark    Input Parameters:
1296c9567895SMark +  pc - the preconditioner context
1297c9567895SMark .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1298c9567895SMark -  n - number of values provided in array
1299c9567895SMark 
1300c9567895SMark    Options Database Key:
1301147403d9SBarry Smith .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1302c9567895SMark 
1303c9567895SMark    Level: intermediate
1304c9567895SMark 
1305c9567895SMark .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1306c9567895SMark @*/
1307c9567895SMark PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1308c9567895SMark {
1309c9567895SMark   PetscFunctionBegin;
1310c9567895SMark   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1311c9567895SMark   if (n) PetscValidIntPointer(v,2);
13129566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n)));
1313c9567895SMark   PetscFunctionReturn(0);
1314c9567895SMark }
1315c9567895SMark 
1316c9567895SMark static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1317c9567895SMark {
1318c9567895SMark   PC_MG   *mg      = (PC_MG*)pc->data;
1319c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1320c9567895SMark   PetscInt i;
1321c9567895SMark   PetscFunctionBegin;
1322c9567895SMark   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1323c9567895SMark   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1324c9567895SMark   PetscFunctionReturn(0);
1325c9567895SMark }
1326c9567895SMark 
1327c9567895SMark /*@
1328c1eae691SMark Adams    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1329c1eae691SMark Adams 
1330c1eae691SMark Adams    Not collective on PC
1331c1eae691SMark Adams 
1332c1eae691SMark Adams    Input Parameters:
1333c1eae691SMark Adams +  pc - the preconditioner context
1334147403d9SBarry Smith -  scale - the threshold value reduction, usually < 1.0
1335c1eae691SMark Adams 
1336c1eae691SMark Adams    Options Database Key:
1337147403d9SBarry Smith .  -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1338c1eae691SMark Adams 
1339055c8bd0SJed Brown    Notes:
1340055c8bd0SJed Brown    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1341055c8bd0SJed Brown    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1342055c8bd0SJed Brown 
1343c1eae691SMark Adams    Level: advanced
1344c1eae691SMark Adams 
1345055c8bd0SJed Brown .seealso: PCGAMGSetThreshold()
1346c1eae691SMark Adams @*/
1347c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1348c1eae691SMark Adams {
13493542efc5SMark F. Adams   PetscFunctionBegin;
1350c1eae691SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13519566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v)));
1352c1eae691SMark Adams   PetscFunctionReturn(0);
1353c1eae691SMark Adams }
1354c1eae691SMark Adams 
1355c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1356c1eae691SMark Adams {
1357c1eae691SMark Adams   PC_MG   *mg      = (PC_MG*)pc->data;
1358c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1359c1eae691SMark Adams   PetscFunctionBegin;
1360c1eae691SMark Adams   pc_gamg->threshold_scale = v;
13613542efc5SMark F. Adams   PetscFunctionReturn(0);
13623542efc5SMark F. Adams }
13633542efc5SMark F. Adams 
1364e20c40e8SBarry Smith /*@C
1365c60c7ad4SBarry Smith    PCGAMGSetType - Set solution method
1366676e1743SMark F. Adams 
1367676e1743SMark F. Adams    Collective on PC
1368676e1743SMark F. Adams 
1369676e1743SMark F. Adams    Input Parameters:
1370c60c7ad4SBarry Smith +  pc - the preconditioner context
1371c60c7ad4SBarry Smith -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1372676e1743SMark F. Adams 
1373676e1743SMark F. Adams    Options Database Key:
1374cab9ed1eSBarry Smith .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1375676e1743SMark F. Adams 
1376676e1743SMark F. Adams    Level: intermediate
1377676e1743SMark F. Adams 
1378cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1379676e1743SMark F. Adams @*/
138019fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1381676e1743SMark F. Adams {
1382676e1743SMark F. Adams   PetscFunctionBegin;
1383676e1743SMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13849566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type)));
1385676e1743SMark F. Adams   PetscFunctionReturn(0);
1386676e1743SMark F. Adams }
1387676e1743SMark F. Adams 
1388e20c40e8SBarry Smith /*@C
1389c60c7ad4SBarry Smith    PCGAMGGetType - Get solution method
1390c60c7ad4SBarry Smith 
1391c60c7ad4SBarry Smith    Collective on PC
1392c60c7ad4SBarry Smith 
1393c60c7ad4SBarry Smith    Input Parameter:
1394c60c7ad4SBarry Smith .  pc - the preconditioner context
1395c60c7ad4SBarry Smith 
1396c60c7ad4SBarry Smith    Output Parameter:
1397c60c7ad4SBarry Smith .  type - the type of algorithm used
1398c60c7ad4SBarry Smith 
1399c60c7ad4SBarry Smith    Level: intermediate
1400c60c7ad4SBarry Smith 
14011c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType
1402c60c7ad4SBarry Smith @*/
1403c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1404c60c7ad4SBarry Smith {
1405c60c7ad4SBarry Smith   PetscFunctionBegin;
1406c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14079566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type)));
1408c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1409c60c7ad4SBarry Smith }
1410c60c7ad4SBarry Smith 
1411c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1412c60c7ad4SBarry Smith {
1413c60c7ad4SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
1414c60c7ad4SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1415c60c7ad4SBarry Smith 
1416c60c7ad4SBarry Smith   PetscFunctionBegin;
1417c60c7ad4SBarry Smith   *type = pc_gamg->type;
1418c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1419c60c7ad4SBarry Smith }
1420c60c7ad4SBarry Smith 
14211e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1422676e1743SMark F. Adams {
14231ab5ffc9SJed Brown   PC_MG          *mg      = (PC_MG*)pc->data;
14241ab5ffc9SJed Brown   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
14255f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PC);
1426676e1743SMark F. Adams 
1427676e1743SMark F. Adams   PetscFunctionBegin;
1428c60c7ad4SBarry Smith   pc_gamg->type = type;
14299566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(GAMGList,type,&r));
14305f80ce2aSJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
14311ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
14329566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
14339566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps)));
1434e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
14353ae0bb68SMark Adams     /* cleaning up common data in pc_gamg - this should disapear someday */
14363ae0bb68SMark Adams     pc_gamg->data_cell_cols = 0;
14373ae0bb68SMark Adams     pc_gamg->data_cell_rows = 0;
14383ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
14393ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
14409566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
14413ae0bb68SMark Adams     pc_gamg->data_sz = 0;
14421ab5ffc9SJed Brown   }
14439566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
14449566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name));
14459566063dSJacob Faibussowitsch   PetscCall((*r)(pc));
1446676e1743SMark F. Adams   PetscFunctionReturn(0);
1447676e1743SMark F. Adams }
1448676e1743SMark F. Adams 
14495adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
14505adeb434SBarry Smith {
14515adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
14525adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1453e7d4b4cbSMark Adams   PetscReal       gc=0, oc=0;
145490db8557SMark Adams 
14555adeb434SBarry Smith   PetscFunctionBegin;
14569566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n"));
14579566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level ="));
14589566063dSJacob Faibussowitsch   for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]));
14599566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
14609566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale));
1461cab9ed1eSBarry Smith   if (pc_gamg->use_aggs_in_asm) {
14629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n"));
1463cab9ed1eSBarry Smith   }
1464171cca9aSMark Adams   if (pc_gamg->use_parallel_coarse_grid_solver) {
14659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n"));
1466171cca9aSMark Adams   }
1467ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1468ce7c7f2fSMark Adams   if (pc_gamg->cpu_pin_coarse_grids) {
14699566063dSJacob Faibussowitsch     /* PetscCall(PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n")); */
1470ce7c7f2fSMark Adams   }
1471ce7c7f2fSMark Adams #endif
1472ce7c7f2fSMark Adams   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
14739566063dSJacob Faibussowitsch   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */
1474ce7c7f2fSMark Adams   /* } else { */
14759566063dSJacob Faibussowitsch   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */
1476ce7c7f2fSMark Adams   /* } */
14775adeb434SBarry Smith   if (pc_gamg->ops->view) {
14789566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->view)(pc,viewer));
14795adeb434SBarry Smith   }
14809566063dSJacob Faibussowitsch   PetscCall(PCMGGetGridComplexity(pc,&gc,&oc));
14819566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g    operator = %g\n",gc,oc));
14825adeb434SBarry Smith   PetscFunctionReturn(0);
14835adeb434SBarry Smith }
14845adeb434SBarry Smith 
14854416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
14865b89ad90SMark F. Adams {
1487676e1743SMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1488676e1743SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
14897e6512fdSJed Brown   PetscBool      flag;
14903b4367a7SBarry Smith   MPI_Comm       comm;
149118c3aa7eSMark   char           prefix[256],tname[32];
1492c1eae691SMark Adams   PetscInt       i,n;
149314a9496bSBarry Smith   const char     *pcpre;
14940a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
14955b89ad90SMark F. Adams   PetscFunctionBegin;
14969566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
14979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"GAMG options"));
14989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
1499bd94a7aaSJed Brown   if (flag) {
15009566063dSJacob Faibussowitsch     PetscCall(PCGAMGSetType(pc,tname));
15011ab5ffc9SJed Brown   }
15029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL));
15039566063dSJacob Faibussowitsch   PetscCall(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));
15049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL));
15059566063dSJacob Faibussowitsch   PetscCall(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));
15069566063dSJacob Faibussowitsch   PetscCall(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));
15079566063dSJacob Faibussowitsch   PetscCall(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));
15089566063dSJacob Faibussowitsch   PetscCall(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));
15099566063dSJacob Faibussowitsch   PetscCall(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));
15109566063dSJacob Faibussowitsch   PetscCall(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));
15119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL));
151218c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
15139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag));
151418c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1515efd3c5ceSMark Adams     if (!flag) n = 1;
1516c1eae691SMark Adams     i = n;
151718c3aa7eSMark     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1518c1eae691SMark Adams   }
1519c9567895SMark   n = PETSC_MG_MAXLEVELS;
15209566063dSJacob Faibussowitsch   PetscCall(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));
1521c9567895SMark   if (!flag) i = 0;
1522c9567895SMark   else i = n;
1523c9567895SMark   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
15249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL));
152518c3aa7eSMark   {
152618c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
152718c3aa7eSMark     n = 2;
15289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag));
152918c3aa7eSMark     if (flag) {
15302c71b3e2SJacob Faibussowitsch       PetscCheckFalse(n != 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
15319566063dSJacob Faibussowitsch       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
153218c3aa7eSMark     }
153318c3aa7eSMark   }
1534b7cbab4eSMark Adams   /* set options for subtype */
15359566063dSJacob Faibussowitsch   if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc));
153618c3aa7eSMark 
15379566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
15389566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
15399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
15405b89ad90SMark F. Adams   PetscFunctionReturn(0);
15415b89ad90SMark F. Adams }
15425b89ad90SMark F. Adams 
15435b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
15445b89ad90SMark F. Adams /*MC
15451cc46a46SBarry Smith      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
15465b89ad90SMark F. Adams 
1547280d9858SJed Brown    Options Database Keys:
1548cab9ed1eSBarry Smith +   -pc_gamg_type <type> - one of agg, geo, or classical
1549cab9ed1eSBarry Smith .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1550cab9ed1eSBarry Smith .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1551cab9ed1eSBarry 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
1552cab9ed1eSBarry 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>
1553cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
1554cab9ed1eSBarry Smith .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
15556008e27bSRichard Tran Mills .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1556c1eae691SMark Adams -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1557cab9ed1eSBarry Smith 
1558cab9ed1eSBarry Smith    Options Database Keys for default Aggregation:
1559cab9ed1eSBarry Smith +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1560cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1561cab9ed1eSBarry Smith -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1562cab9ed1eSBarry Smith 
1563db9745e2SBarry Smith    Multigrid options:
1564db9745e2SBarry Smith +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1565db9745e2SBarry Smith .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1566db9745e2SBarry Smith .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1567cab9ed1eSBarry Smith -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
15685b89ad90SMark F. Adams 
156995452b02SPatrick Sanan   Notes:
157095452b02SPatrick Sanan     In order to obtain good performance for PCGAMG for vector valued problems you must
1571db9745e2SBarry Smith        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1572db9745e2SBarry Smith        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1573db9745e2SBarry Smith        See the Users Manual Chapter 4 for more details
15741cc46a46SBarry Smith 
15755b89ad90SMark F. Adams   Level: intermediate
1576280d9858SJed Brown 
15771cc46a46SBarry Smith .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
157873f7197eSJed Brown            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig()
15795b89ad90SMark F. Adams M*/
1580b2573a8aSBarry Smith 
15818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
15825b89ad90SMark F. Adams {
15835b89ad90SMark F. Adams   PC_GAMG *pc_gamg;
15845b89ad90SMark F. Adams   PC_MG   *mg;
15855b89ad90SMark F. Adams 
15865b89ad90SMark F. Adams   PetscFunctionBegin;
15871c1aac46SBarry Smith    /* register AMG type */
15889566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
15891c1aac46SBarry Smith 
15905b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
15919566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
15929566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
15935b89ad90SMark F. Adams 
15945b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
15959566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg));
15969566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
15975b89ad90SMark F. Adams   mg           = (PC_MG*)pc->data;
15985b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
15995b89ad90SMark F. Adams 
16009566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg->ops));
16011ab5ffc9SJed Brown 
16029d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
16039d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
16040a545947SLisandro Dalcin   pc_gamg->data    = NULL;
16055b89ad90SMark F. Adams 
16069d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
16075b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
16085b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
16095b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
16105b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
16115adeb434SBarry Smith   mg->view                = PCView_GAMG;
16125b89ad90SMark F. Adams 
16139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG));
16149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG));
16159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG));
16169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG));
16179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG));
16189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG));
16199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG));
16209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG));
16219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG));
16229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG));
16239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG));
16249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG));
16259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG));
16269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG));
16279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG));
16289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG));
16299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG));
16309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG));
16319d5b6da9SMark F. Adams   pc_gamg->repart           = PETSC_FALSE;
1632d3042614SMark Adams   pc_gamg->reuse_prol       = PETSC_FALSE;
16330c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1634171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1635a0095786SMark   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1636a0095786SMark   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1637038f3aa4SMark F. Adams   pc_gamg->min_eq_proc      = 50;
163825a145a7SMark Adams   pc_gamg->coarse_eq_limit  = 50;
16399566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS));
1640c1eae691SMark Adams   pc_gamg->threshold_scale = 1.;
164118c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
16429ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
16437e6512fdSJed Brown   pc_gamg->use_sa_esteig    = PETSC_TRUE;
164418c3aa7eSMark   pc_gamg->emin             = 0;
164518c3aa7eSMark   pc_gamg->emax             = 0;
164618c3aa7eSMark 
1647c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
16489d5b6da9SMark F. Adams 
1649bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
16509566063dSJacob Faibussowitsch   PetscCall(PCGAMGSetType(pc,PCGAMGAGG));
16515b89ad90SMark F. Adams   PetscFunctionReturn(0);
16525b89ad90SMark F. Adams }
16533e3471ccSMark Adams 
16543e3471ccSMark Adams /*@C
16553e3471ccSMark Adams  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
16568a690491SBarry Smith     from PCInitializePackage().
16573e3471ccSMark Adams 
16583e3471ccSMark Adams  Level: developer
16593e3471ccSMark Adams 
16603e3471ccSMark Adams  .seealso: PetscInitialize()
16613e3471ccSMark Adams @*/
16623e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void)
16633e3471ccSMark Adams {
16644555aa8cSStefano Zampini   PetscInt       l;
16653e3471ccSMark Adams 
16663e3471ccSMark Adams   PetscFunctionBegin;
16673e3471ccSMark Adams   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
16683e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
16699566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO));
16709566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG));
16719566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical));
16729566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1673c1c463dbSMark Adams 
1674c1c463dbSMark Adams   /* general events */
16759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG));
16769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO));
16779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG));
16789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO));
16799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG));
16809566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO));
16819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG));
1682c1c463dbSMark Adams 
16839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]));
16849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  Create Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]));
16859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  Filter Graph", PC_CLASSID, &petsc_gamg_setup_events[SET16]));
16865b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
16875b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
16885b89ad90SMark F. Adams   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
16899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]));
16909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]));
16919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]));
16929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]));
16939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]));
16949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]));
16959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]));
16969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]));
16979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]));
16989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]));
16999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]));
17009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]));
17014555aa8cSStefano Zampini   for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
17024555aa8cSStefano Zampini     char ename[32];
17035b89ad90SMark F. Adams 
17049566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l));
17059566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
17069566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l));
17079566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
17089566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l));
17099566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
17104555aa8cSStefano Zampini   }
17115b89ad90SMark F. Adams   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
17125b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
17135b89ad90SMark F. Adams   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
17145b89ad90SMark F. Adams   /* create timer stages */
17154555aa8cSStefano Zampini #if defined(GAMG_STAGES)
17165b89ad90SMark F. Adams   {
17175b89ad90SMark F. Adams     char     str[32];
17185b89ad90SMark F. Adams     PetscInt lidx;
17195b89ad90SMark F. Adams     sprintf(str,"MG Level %d (finest)",0);
17209566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
17215b89ad90SMark F. Adams     for (lidx=1; lidx<9; lidx++) {
1722c9567895SMark       sprintf(str,"MG Level %d",(int)lidx);
17239566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister(str, &gamg_stages[lidx]));
17245b89ad90SMark F. Adams     }
17255b89ad90SMark F. Adams   }
17265b89ad90SMark F. Adams #endif
17273e3471ccSMark Adams   PetscFunctionReturn(0);
17283e3471ccSMark Adams }
17293e3471ccSMark Adams 
17303e3471ccSMark Adams /*@C
17311c1aac46SBarry Smith  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
17321c1aac46SBarry Smith     called from PetscFinalize() automatically.
17333e3471ccSMark Adams 
17343e3471ccSMark Adams  Level: developer
17353e3471ccSMark Adams 
17363e3471ccSMark Adams  .seealso: PetscFinalize()
17373e3471ccSMark Adams @*/
17383e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void)
17393e3471ccSMark Adams {
17403e3471ccSMark Adams   PetscFunctionBegin;
17413e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
17429566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&GAMGList));
17433e3471ccSMark Adams   PetscFunctionReturn(0);
17443e3471ccSMark Adams }
1745a36cf38bSToby Isaac 
1746a36cf38bSToby Isaac /*@C
1747a36cf38bSToby Isaac  PCGAMGRegister - Register a PCGAMG implementation.
1748a36cf38bSToby Isaac 
1749a36cf38bSToby Isaac  Input Parameters:
1750a36cf38bSToby Isaac  + type - string that will be used as the name of the GAMG type.
1751a36cf38bSToby Isaac  - create - function for creating the gamg context.
1752a36cf38bSToby Isaac 
1753a36cf38bSToby Isaac   Level: advanced
1754a36cf38bSToby Isaac 
17551c1aac46SBarry Smith  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1756a36cf38bSToby Isaac @*/
1757a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1758a36cf38bSToby Isaac {
1759a36cf38bSToby Isaac   PetscFunctionBegin;
17609566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
17619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList,type,create));
1762a36cf38bSToby Isaac   PetscFunctionReturn(0);
1763a36cf38bSToby Isaac }
1764