xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
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 
16849bee69SMark Adams PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET];
174555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
180cbbd2e1SMark F. Adams 
19849bee69SMark Adams // #define GAMG_STAGES
204555aa8cSStefano Zampini #if defined(GAMG_STAGES)
2118c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
22b4fbaa2aSMark F. Adams #endif
23f96513f1SMatthew G Knepley 
240a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
253e3471ccSMark Adams static PetscBool         PCGAMGPackageInitialized;
269d5b6da9SMark F. Adams 
2766976f2fSJacob Faibussowitsch static PetscErrorCode PCReset_GAMG(PC pc)
28d71ae5a4SJacob Faibussowitsch {
29d3d6bff4SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
30d3d6bff4SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
31d3d6bff4SMark F. Adams 
32d3d6bff4SMark F. Adams   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
341c1aac46SBarry Smith   pc_gamg->data_sz = 0;
359566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->orig_data));
365f80ce2aSJacob Faibussowitsch   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) {
3718c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
3818c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
3918c3aa7eSMark   }
4018c3aa7eSMark   pc_gamg->emin = 0;
4118c3aa7eSMark   pc_gamg->emax = 0;
42978e3cbaSStefano Zampini   PetscCall(PCReset_MG(pc));
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44a2f3521dSMark F. Adams }
45a2f3521dSMark F. Adams 
465b89ad90SMark F. Adams /*
47c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
48a147abb0SMark F. Adams      of active processors.
495b89ad90SMark F. Adams 
505b89ad90SMark F. Adams    Input Parameter:
51a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
52a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
539d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
54c5bfad50SMark F. Adams    . cr_bs - coarse block size
553530afc2SMark F. Adams    In/Output Parameter:
56a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
57afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5811e60469SMark F. Adams    Output Parameter:
593530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
605b89ad90SMark F. Adams */
61d71ae5a4SJacob Faibussowitsch 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)
62d71ae5a4SJacob Faibussowitsch {
639d5b6da9SMark F. Adams   PC_MG      *mg      = (PC_MG *)pc->data;
64486a8d0bSJed Brown   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
65a2f3521dSMark F. Adams   Mat         Cmat, Pold = *a_P_inout;
663b4367a7SBarry Smith   MPI_Comm    comm;
67c5df96a5SBarry Smith   PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc;
683ae0bb68SMark Adams   PetscInt    ncrs_eq, ncrs, f_bs;
695b89ad90SMark F. Adams 
705b89ad90SMark F. Adams   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm));
729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
749566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
75849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
779566063dSJacob Faibussowitsch   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
79849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
80038e3b61SMark F. Adams 
81ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
82ce7c7f2fSMark Adams 
833ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
849566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
853ae0bb68SMark Adams   if (pc_gamg->data_cell_rows > 0) {
863ae0bb68SMark Adams     ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows;
8773911c69SBarry Smith   } else {
883ae0bb68SMark Adams     PetscInt bs;
899566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(Cmat, &bs));
903ae0bb68SMark Adams     ncrs = ncrs_eq / bs;
913ae0bb68SMark Adams   }
92c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
93c9567895SMark   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. */
94c9567895SMark #if defined(PETSC_HAVE_CUDA)
95c9567895SMark     PetscShmComm pshmcomm;
96c9567895SMark     PetscMPIInt  locrank;
97c9567895SMark     MPI_Comm     loccomm;
98c9567895SMark     PetscInt     s_nnodes, r_nnodes, new_new_size;
99c9567895SMark     cudaError_t  cerr;
100c9567895SMark     int          devCount;
1019566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGet(comm, &pshmcomm));
1029566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm));
1039566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
104c9567895SMark     s_nnodes = !locrank;
105712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm));
10663a3b9bcSJacob Faibussowitsch     PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes);
107c9567895SMark     devCount = 0;
108c9567895SMark     cerr     = cudaGetDeviceCount(&devCount);
109c9567895SMark     cudaGetLastError();                         /* Reset the last error */
110c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
111c9567895SMark       new_new_size = r_nnodes * devCount;
112c9567895SMark       new_size     = new_new_size;
11363a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Fine grid with Cuda. %" PetscInt_FMT " nodes. Change new active set size %d --> %d (devCount=%d #nodes=%" PetscInt_FMT ")\n", ((PetscObject)pc)->prefix, r_nnodes, nactive, new_size, devCount, r_nnodes));
114c9567895SMark     } else {
1159d3446b2SPierre Jolivet       PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.\n", ((PetscObject)pc)->prefix));
116c9567895SMark       goto HEURISTIC;
117c9567895SMark     }
118c9567895SMark #else
119c9567895SMark     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here");
120c9567895SMark #endif
121c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
1227509f629SMark Adams     if (nactive < pc_gamg->level_reduction_factors[pc_gamg->current_level]) {
1237509f629SMark Adams       new_size = 1;
1247509f629SMark Adams       PetscCall(PetscInfo(pc, "%s: reduction factor too small for %d active processes: reduce to one process\n", ((PetscObject)pc)->prefix, new_size));
1257509f629SMark Adams     } else {
12663a3b9bcSJacob Faibussowitsch       PetscCheck(nactive % pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of active process %d wrt reduction factor %" PetscInt_FMT, nactive, pc_gamg->level_reduction_factors[pc_gamg->current_level]);
127c9567895SMark       new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level];
12863a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Manually setting reduction to %d active processes (%d/%" PetscInt_FMT ")\n", ((PetscObject)pc)->prefix, new_size, nactive, pc_gamg->level_reduction_factors[pc_gamg->current_level]));
1297509f629SMark Adams     }
130c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
131c9567895SMark     new_size = 1;
1329566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size));
133c9567895SMark   } else {
134472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
135c9567895SMark #if defined(PETSC_HAVE_CUDA)
136c9567895SMark   HEURISTIC:
137c9567895SMark #endif
1389566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
139a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
140da81f932SPierre Jolivet     if (!new_size) new_size = 1;                                                       /* not likely, possible? */
141c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive;                                  /* no change, rare */
1429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size));
143a2f3521dSMark F. Adams   }
1442e3501ffSMark Adams   if (new_size == nactive) {
145ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
146ce7c7f2fSMark Adams     if (new_size < size) {
147ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
1489566063dSJacob 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));
149ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
1509566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE));
1519566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE));
152ce7c7f2fSMark Adams       }
153ce7c7f2fSMark Adams     }
154ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1552e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
156192c0e8bSMark Adams     PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1;
157885364a3SMark Adams     IS        is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices;
158849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
15971959b99SBarry Smith     nloc_old = ncrs_eq / cr_bs;
16063a3b9bcSJacob Faibussowitsch     PetscCheck(ncrs_eq % cr_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ncrs_eq %" PetscInt_FMT " not divisible by cr_bs %" PetscInt_FMT, ncrs_eq, cr_bs);
161ce7c7f2fSMark Adams     /* get new_size and rfactor */
162ce7c7f2fSMark Adams     if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
163ce7c7f2fSMark Adams       /* find factor */
164ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
165ce7c7f2fSMark Adams       else {
166ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
167ce7c7f2fSMark Adams         jj                  = -1;
168ce7c7f2fSMark Adams         for (kk = 1; kk <= size; kk++) {
169ce7c7f2fSMark Adams           if (!(size % kk)) { /* a candidate */
170ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size;
171ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */
172ce7c7f2fSMark Adams             if (fact > best_fact) {
1739371c9d4SSatish Balay               best_fact = fact;
1749371c9d4SSatish Balay               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;
18663a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, new_size, ncrs_eq));
187849bee69SMark Adams         PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
1883ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
189ce7c7f2fSMark Adams       }
190ce7c7f2fSMark Adams     }
191a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
1929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &counts));
193849bee69SMark Adams     if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
1945a9b9e01SMark F. Adams       Mat adj;
195849bee69SMark Adams       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
19663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Repartition: size (active): %d --> %d, %" PetscInt_FMT " 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) {
2449371c9d4SSatish Balay           PetscViewer viewer;
2459371c9d4SSatish Balay           char        fname[32];
24663a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev));
2473ba16761SJacob Faibussowitsch           PetscCall(PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer));
2489566063dSJacob Faibussowitsch           PetscCall(MatView(tMat, viewer));
2499566063dSJacob Faibussowitsch           PetscCall(PetscViewerDestroy(&viewer));
250b4fbaa2aSMark F. Adams         }
2519566063dSJacob Faibussowitsch         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
2529566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tMat));
253a2f3521dSMark F. Adams       } /* create 'adj' */
254f150b916SMark F. Adams 
255a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2565a9b9e01SMark F. Adams         char            prefix[256];
2575a9b9e01SMark F. Adams         const char     *pcpre;
258b4fbaa2aSMark F. Adams         const PetscInt *is_idx;
259b4fbaa2aSMark F. Adams         MatPartitioning mpart;
260a4b7d37bSMark F. Adams         IS              proc_is;
2612f03bc48SMark F. Adams 
2629566063dSJacob Faibussowitsch         PetscCall(MatPartitioningCreate(comm, &mpart));
2639566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
2649566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
2659566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
2669566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
2679566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetFromOptions(mpart));
2689566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, new_size));
2699566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &proc_is));
2709566063dSJacob Faibussowitsch         PetscCall(MatPartitioningDestroy(&mpart));
2715a9b9e01SMark F. Adams 
2725ef31b24SMark F. Adams         /* collect IS info */
2739566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
2749566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(proc_is, &is_idx));
275a2f3521dSMark F. Adams         for (kk = jj = 0; kk < nloc_old; kk++) {
2769371c9d4SSatish Balay           for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ }
2775ef31b24SMark F. Adams         }
2789566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(proc_is, &is_idx));
2799566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&proc_is));
2805ef31b24SMark F. Adams       }
2819566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&adj));
2825a9b9e01SMark F. Adams 
2839566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
2849566063dSJacob Faibussowitsch       PetscCall(PetscFree(newproc_idx));
285849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
28631cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
287ce7c7f2fSMark Adams       PetscInt targetPE;
28808401ef6SPierre Jolivet       PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen");
28963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " 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));
3086aad120cSJose E. Roman     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxiliary 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 */
38111e60469SMark F. Adams     /*
3827dae84e0SHong Zhang       Invert for MatCreateSubMatrix
38311e60469SMark F. Adams     */
3849566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
3859566063dSJacob Faibussowitsch     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
3869566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
3879371c9d4SSatish Balay     if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ }
3883cb8563fSToby Isaac     if (Pcolumnperm) {
3899566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
3903cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3913cb8563fSToby Isaac     }
3929566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_num));
393849bee69SMark Adams 
394a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
395a2f3521dSMark F. Adams     {
396a2f3521dSMark F. Adams       Mat       mat;
397b94d7dedSBarry Smith       PetscBool isset, isspd, isher;
39890db8557SMark Adams #if !defined(PETSC_USE_COMPLEX)
399b94d7dedSBarry Smith       PetscBool issym;
400b94d7dedSBarry Smith #endif
401b94d7dedSBarry Smith 
402b94d7dedSBarry Smith       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
403b94d7dedSBarry Smith       PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ?
404b94d7dedSBarry Smith       if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd));
405b94d7dedSBarry Smith       else {
406b94d7dedSBarry Smith         PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher));
407b94d7dedSBarry Smith         if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher));
408b94d7dedSBarry Smith         else {
409b94d7dedSBarry Smith #if !defined(PETSC_USE_COMPLEX)
410b94d7dedSBarry Smith           PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym));
411b94d7dedSBarry Smith           if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym));
41290db8557SMark Adams #endif
41390db8557SMark Adams         }
41490db8557SMark Adams       }
415a2f3521dSMark F. Adams       *a_Amat_crs = mat;
416a2f3521dSMark F. Adams     }
4179566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cmat));
418a2f3521dSMark F. Adams 
41911e60469SMark F. Adams     /* prolongator */
42011e60469SMark F. Adams     {
42111e60469SMark F. Adams       IS       findices;
422a2f3521dSMark F. Adams       PetscInt Istart, Iend;
423a2f3521dSMark F. Adams       Mat      Pnew;
42462294041SBarry Smith 
4259566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
4269566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices));
4279566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(findices, f_bs));
4289566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
4299566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&findices));
4309566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
431c5bfad50SMark F. Adams 
4329566063dSJacob Faibussowitsch       PetscCall(MatDestroy(a_P_inout));
433a2f3521dSMark F. Adams 
434a2f3521dSMark F. Adams       /* output - repartitioned */
435a2f3521dSMark F. Adams       *a_P_inout = Pnew;
436e33ef3b1SMark F. Adams     }
4379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&new_eq_indices));
4385b89ad90SMark F. Adams 
439c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
440ce7c7f2fSMark Adams 
441ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
442ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
443ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4448bca76a6SMark Adams       static PetscInt llev = 2;
44563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++));
446ce7c7f2fSMark Adams #endif
4479566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE));
4489566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE));
449adf5291fSStefano Zampini       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
450ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
451ce7c7f2fSMark Adams         PetscMPIInt size;
4529566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
453ce7c7f2fSMark Adams         if (size > 1) {
454ce7c7f2fSMark Adams           Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data;
4559566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE));
4569566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE));
457ce7c7f2fSMark Adams         }
458ce7c7f2fSMark Adams       }
459ce7c7f2fSMark Adams     }
460849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
461849bee69SMark Adams   } // processor reduce
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4635b89ad90SMark F. Adams }
4645b89ad90SMark F. Adams 
465bae903cbSmarkadams4 // used in GEO
466d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2)
467d71ae5a4SJacob Faibussowitsch {
4684b1575e2SStefano Zampini   const char *prefix;
4694b1575e2SStefano Zampini   char        addp[32];
4704b1575e2SStefano Zampini   PC_MG      *mg      = (PC_MG *)a_pc->data;
4714b1575e2SStefano Zampini   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
4724b1575e2SStefano Zampini 
4734b1575e2SStefano Zampini   PetscFunctionBegin;
4749566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(a_pc, &prefix));
47563a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1));
4769566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2));
4779566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(*Gmat2, prefix));
47863a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level));
4799566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(*Gmat2, addp));
480b94d7dedSBarry Smith   if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) {
4819566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB));
482b4da3a1bSStefano Zampini   } else {
4839566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
4849566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB));
485b4da3a1bSStefano Zampini   }
4869566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(*Gmat2));
4879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0));
4889566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(*Gmat2));
4899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0));
4909566063dSJacob Faibussowitsch   PetscCall(MatProductClear(*Gmat2));
4914b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
4924b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4944b1575e2SStefano Zampini }
4954b1575e2SStefano Zampini 
4965b89ad90SMark F. Adams /*
4975b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4985b89ad90SMark F. Adams                     by setting data structures and options.
4995b89ad90SMark F. Adams 
5005b89ad90SMark F. Adams    Input Parameter:
5015b89ad90SMark F. Adams .  pc - the preconditioner context
5025b89ad90SMark F. Adams 
5035b89ad90SMark F. Adams */
50466976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_GAMG(PC pc)
505d71ae5a4SJacob Faibussowitsch {
5069d5b6da9SMark F. Adams   PC_MG      *mg      = (PC_MG *)pc->data;
5075b89ad90SMark F. Adams   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
5082adcac29SMark F. Adams   Mat         Pmat    = pc->pmat;
50918c3aa7eSMark   PetscInt    fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS];
5103b4367a7SBarry Smith   MPI_Comm    comm;
511c5df96a5SBarry Smith   PetscMPIInt rank, size, nactivepe;
51218c3aa7eSMark   Mat         Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS];
51318c3aa7eSMark   IS         *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
5144279555eSSatish Balay   PetscBool   is_last = PETSC_FALSE;
5154279555eSSatish Balay #if defined(PETSC_USE_INFO)
516a2f3521dSMark F. Adams   PetscLogDouble nnz0 = 0., nnztot = 0.;
517569f4572SMark Adams   MatInfo        info;
5184279555eSSatish Balay #endif
5195ef31b24SMark F. Adams 
5205b89ad90SMark F. Adams   PetscFunctionBegin;
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
5229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
524849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
5258abdc6daSStefano Zampini   if (pc->setupcalled) {
5268abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
527878e152fSMark F. Adams       /* reset everything */
5289566063dSJacob Faibussowitsch       PetscCall(PCReset_MG(pc));
529878e152fSMark F. Adams       pc->setupcalled = 0;
530806fa848SBarry Smith     } else {
53184d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
53203a628feSMark F. Adams       /* just do Galerkin grids */
53358471d46SMark F. Adams       Mat B, dA, dB;
5349d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
5354555aa8cSStefano Zampini         PetscInt gl;
53658471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5379566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB));
53858471d46SMark F. Adams         /* (re)set to get dirty flag */
5399566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB));
54058471d46SMark F. Adams 
5414555aa8cSStefano Zampini         for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) {
5428abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX;
543849bee69SMark Adams #if defined(GAMG_STAGES)
544849bee69SMark Adams           PetscCall(PetscLogStagePush(gamg_stages[gl]));
545849bee69SMark Adams #endif
5468abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
5479566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B));
5488abdc6daSStefano Zampini           if (B->product) {
549ad540459SPierre Jolivet             if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX;
5508abdc6daSStefano Zampini           }
5519566063dSJacob Faibussowitsch           if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A));
5528abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
55363a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
5548abdc6daSStefano Zampini           } else {
55563a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
5568abdc6daSStefano Zampini           }
5579566063dSJacob Faibussowitsch           PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0));
5589566063dSJacob Faibussowitsch           PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B));
5599566063dSJacob Faibussowitsch           PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0));
56063b77682SMark Adams           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
5619566063dSJacob Faibussowitsch           PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B));
562e1cf1444SMark Adams           // check for redoing eigen estimates
563e1cf1444SMark Adams           if (pc_gamg->recompute_esteig) {
564e1cf1444SMark Adams             PetscBool ischeb;
565e1cf1444SMark Adams             KSP       smoother;
566e1cf1444SMark Adams             PetscCall(PCMGGetSmoother(pc, level + 1, &smoother));
567e1cf1444SMark Adams             PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
568e1cf1444SMark Adams             if (ischeb) {
569e1cf1444SMark Adams               KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
570e1cf1444SMark Adams               cheb->emin_provided = 0;
571e1cf1444SMark Adams               cheb->emax_provided = 0;
572e1cf1444SMark Adams             }
573e1cf1444SMark Adams             /* we could call PetscCall(KSPChebyshevSetEigenvalues(smoother, 0, 0)); but the logic does not work properly */
574e1cf1444SMark Adams           }
575e1cf1444SMark Adams           // inc
57658471d46SMark F. Adams           dB = B;
577849bee69SMark Adams #if defined(GAMG_STAGES)
578849bee69SMark Adams           PetscCall(PetscLogStagePop());
579849bee69SMark Adams #endif
58058471d46SMark F. Adams         }
5815f8cf99dSMark F. Adams       }
582d5280255SMark F. Adams 
5839566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
584849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
5853ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
586eb07cef2SMark F. Adams     }
587878e152fSMark F. Adams   }
588f6536408SMark F. Adams 
589878e152fSMark F. Adams   if (!pc_gamg->data) {
590878e152fSMark F. Adams     if (pc_gamg->orig_data) {
5919566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(Pmat, &bs));
5929566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
5932fa5cd67SKarl Rupp 
594878e152fSMark F. Adams       pc_gamg->data_sz        = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols;
595878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
596878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5972fa5cd67SKarl Rupp 
5989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
599878e152fSMark F. Adams       for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
600806fa848SBarry Smith     } else {
6015f80ce2aSJacob Faibussowitsch       PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data");
6029566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat));
6039d5b6da9SMark F. Adams     }
604878e152fSMark F. Adams   }
605878e152fSMark F. Adams 
606878e152fSMark F. Adams   /* cache original data for reuse */
6071c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
6089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
609878e152fSMark F. Adams     for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
610878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
611878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
612878e152fSMark F. Adams   }
613038e3b61SMark F. Adams 
614302f38e8SMark F. Adams   /* get basic dims */
6159566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Pmat, &bs));
6169566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Pmat, &M, &N));
61784d3f75bSMark F. Adams 
6184279555eSSatish Balay #if defined(PETSC_USE_INFO)
6199566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */
620569f4572SMark Adams   nnz0   = info.nz_used;
621569f4572SMark Adams   nnztot = info.nz_used;
6224279555eSSatish Balay #endif
623bae903cbSmarkadams4   PetscCall(PetscInfo(pc, "%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", np=%d\n", ((PetscObject)pc)->prefix, 0, M, pc_gamg->data_cell_rows, pc_gamg->data_cell_cols, (PetscInt)(nnz0 / (PetscReal)M + 0.5), size));
624569f4572SMark Adams 
625a2f3521dSMark F. Adams   /* Get A_i and R_i */
62662294041SBarry Smith   for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) {
6279ab59c8bSMark Adams     pc_gamg->current_level = level;
62863a3b9bcSJacob Faibussowitsch     PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level);
6295b89ad90SMark F. Adams     level1 = level + 1;
6304555aa8cSStefano Zampini #if defined(GAMG_STAGES)
631849bee69SMark Adams     if (!gamg_stages[level]) {
632849bee69SMark Adams       char str[32];
633a364092eSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level));
634849bee69SMark Adams       PetscCall(PetscLogStageRegister(str, &gamg_stages[level]));
635849bee69SMark Adams     }
6369566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(gamg_stages[level]));
637b4fbaa2aSMark F. Adams #endif
638c8b0795cSMark F. Adams     { /* construct prolongator */
639725b86d8SJed Brown       Mat               Gmat;
6400cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6417700e67bSMark Adams       Mat               Prol11;
642c8b0795cSMark F. Adams 
6432d776b49SBarry Smith       PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat));
6449566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists));
6459566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11));
646c8b0795cSMark F. Adams 
647a2f3521dSMark F. Adams       /* could have failed to create new level */
648a2f3521dSMark F. Adams       if (Prol11) {
649f7df55f0SStefano Zampini         const char *prefix;
650f7df55f0SStefano Zampini         char        addp[32];
651f7df55f0SStefano Zampini 
6529d5b6da9SMark F. Adams         /* get new block size of coarse matrices */
6539566063dSJacob Faibussowitsch         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs));
654a2f3521dSMark F. Adams 
655fd1112cbSBarry Smith         if (pc_gamg->ops->optprolongator) {
656c8b0795cSMark F. Adams           /* smooth */
6579566063dSJacob Faibussowitsch           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
658c8b0795cSMark F. Adams         }
659c8b0795cSMark F. Adams 
6600c3bc534SBarry Smith         if (pc_gamg->use_aggs_in_asm) {
6611b18a24aSMark Adams           PetscInt bs;
662849bee69SMark Adams           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method
6639566063dSJacob Faibussowitsch           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
664ffc955d6SMark F. Adams         }
665ffc955d6SMark F. Adams 
6669566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
6679566063dSJacob Faibussowitsch         PetscCall(MatSetOptionsPrefix(Prol11, prefix));
6689566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level));
6699566063dSJacob Faibussowitsch         PetscCall(MatAppendOptionsPrefix(Prol11, addp));
67091f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
671f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
6729566063dSJacob Faibussowitsch         PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
6739566063dSJacob Faibussowitsch         PetscCall(MatSetFromOptions(Prol11));
6744bde40a0SMark Adams         Parr[level1] = Prol11;
6754bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
6764bde40a0SMark Adams 
6779566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat));
6789566063dSJacob Faibussowitsch       PetscCall(PetscCDDestroy(agg_lists));
679a2f3521dSMark F. Adams     }                           /* construct prolongator scope */
6807f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
681171cca9aSMark Adams     if (!Parr[level1]) {        /* failed to coarsen */
68263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
6834555aa8cSStefano Zampini #if defined(GAMG_STAGES)
6849566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePop());
685a90e85d9SMark Adams #endif
686c8b0795cSMark F. Adams       break;
687c8b0795cSMark F. Adams     }
6889566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
6892472a847SBarry Smith     PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?");
690171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
6910e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE;
692849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
6939566063dSJacob Faibussowitsch     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
694849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
695a2f3521dSMark F. Adams 
6969566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
6974279555eSSatish Balay #if defined(PETSC_USE_INFO)
6989566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
699569f4572SMark Adams     nnztot += info.nz_used;
7004279555eSSatish Balay #endif
701bae903cbSmarkadams4     PetscCall(PetscInfo(pc, "%s: %d) N=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", %d active pes\n", ((PetscObject)pc)->prefix, (int)level1, M, pc_gamg->data_cell_cols, (PetscInt)(info.nz_used / (PetscReal)M), nactivepe));
702569f4572SMark Adams 
7034555aa8cSStefano Zampini #if defined(GAMG_STAGES)
7049566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
705b4fbaa2aSMark F. Adams #endif
706a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
7079ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) {
70863a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: HARD stop of coarsening on level %" PetscInt_FMT ".  Grid too small: %" PetscInt_FMT " block nodes\n", ((PetscObject)pc)->prefix, level, M / bs));
709a90e85d9SMark Adams       level++;
710a90e85d9SMark Adams       break;
711a90e85d9SMark Adams     }
712c8b0795cSMark F. Adams   } /* levels */
7139566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
714c8b0795cSMark F. Adams 
715ba4ceb06Smarkadams4   PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0));
7169d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7175b89ad90SMark F. Adams   fine_level       = level;
7189566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL));
7195b89ad90SMark F. Adams 
72062294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7210ed2132dSStefano Zampini 
722d5280255SMark F. Adams     /* set default smoothers & set operators */
72362294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) {
724ffc955d6SMark F. Adams       KSP smoother;
725ffc955d6SMark F. Adams       PC  subpc;
726a2f3521dSMark F. Adams 
7279566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7289566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(smoother, &subpc));
729ffc955d6SMark F. Adams 
7309566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
731a2f3521dSMark F. Adams       /* set ops */
7329566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
7339566063dSJacob Faibussowitsch       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1]));
734a2f3521dSMark F. Adams 
735a2f3521dSMark F. Adams       /* set defaults */
7369566063dSJacob Faibussowitsch       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
737a2f3521dSMark F. Adams 
7380c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
7390c3bc534SBarry Smith       if (pc_gamg->use_aggs_in_asm) {
7402d3561bbSSatish Balay         PetscInt sz;
7417a28f3e5SMark Adams         IS      *iss;
742a2f3521dSMark F. Adams 
7432d3561bbSSatish Balay         sz  = nASMBlocksArr[level];
7447a28f3e5SMark Adams         iss = ASMLocalIDsArr[level];
7459566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCASM));
7469566063dSJacob Faibussowitsch         PetscCall(PCASMSetOverlap(subpc, 0));
7479566063dSJacob Faibussowitsch         PetscCall(PCASMSetType(subpc, PC_ASM_BASIC));
7487f66b68fSMark Adams         if (!sz) {
749ffc955d6SMark F. Adams           IS is;
7509566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
7519566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
7529566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
753806fa848SBarry Smith         } else {
754a94c3b12SMark F. Adams           PetscInt kk;
7552eab5cd7Smarkadams4           PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL));
75648a46eb9SPierre Jolivet           for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk]));
7579566063dSJacob Faibussowitsch           PetscCall(PetscFree(iss));
758ffc955d6SMark F. Adams         }
7590298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
760ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
761806fa848SBarry Smith       } else {
7629566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCJACOBI));
763ffc955d6SMark F. Adams       }
764d5280255SMark F. Adams     }
765d5280255SMark F. Adams     {
766d5280255SMark F. Adams       /* coarse grid */
7679371c9d4SSatish Balay       KSP      smoother, *k2;
7689371c9d4SSatish Balay       PC       subpc, pc2;
7699371c9d4SSatish Balay       PetscInt ii, first;
7709371c9d4SSatish Balay       Mat      Lmat = Aarr[(level = pc_gamg->Nlevels - 1)];
7719371c9d4SSatish Balay       lidx          = 0;
7720ed2132dSStefano Zampini 
7739566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7749566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
775cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
7769566063dSJacob Faibussowitsch         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
7779566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(smoother, &subpc));
7789566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCBJACOBI));
7799566063dSJacob Faibussowitsch         PetscCall(PCSetUp(subpc));
7809566063dSJacob Faibussowitsch         PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2));
78163a3b9bcSJacob Faibussowitsch         PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii);
7829566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(k2[0], &pc2));
7839566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc2, PCLU));
7849566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS));
7859566063dSJacob Faibussowitsch         PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1));
7869566063dSJacob Faibussowitsch         PetscCall(KSPSetType(k2[0], KSPPREONLY));
787d5280255SMark F. Adams       }
788cf8ae1d3SMark Adams     }
789d5280255SMark F. Adams 
790d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
791d0609cedSBarry Smith     PetscObjectOptionsBegin((PetscObject)pc);
792dbbe0bcdSBarry Smith     PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject));
793d0609cedSBarry Smith     PetscOptionsEnd();
7949566063dSJacob Faibussowitsch     PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
795d5280255SMark F. Adams 
796f1580f4eSBarry Smith     /* set cheby eigen estimates from SA to use in the solver */
7977e6512fdSJed Brown     if (pc_gamg->use_sa_esteig) {
79818c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) {
79918c3aa7eSMark         KSP       smoother;
80018c3aa7eSMark         PetscBool ischeb;
8010ed2132dSStefano Zampini 
8029566063dSJacob Faibussowitsch         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
8039566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
80418c3aa7eSMark         if (ischeb) {
80518c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
8060ed2132dSStefano Zampini 
8072de708cbSJed Brown           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
808efe053fcSJed Brown           if (mg->max_eigen_DinvA[level] > 0) {
8097e6512fdSJed Brown             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
8107e6512fdSJed Brown             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
81118c3aa7eSMark             PetscReal emax, emin;
8120ed2132dSStefano Zampini 
81318c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
81418c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
81563a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pc, "%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %" PetscInt_FMT " (N=%" PetscInt_FMT ") with emax = %g emin = %g\n", ((PetscObject)pc)->prefix, level, Aarr[level]->rmap->N, (double)emax, (double)emin));
8160a94a983SJed Brown             cheb->emin_provided = emin;
8170a94a983SJed Brown             cheb->emax_provided = emax;
81818c3aa7eSMark           }
81918c3aa7eSMark         }
82018c3aa7eSMark       }
82118c3aa7eSMark     }
8220ed2132dSStefano Zampini 
8239566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8240ed2132dSStefano Zampini 
825d5280255SMark F. Adams     /* clean up */
826d5280255SMark F. Adams     for (level = 1; level < pc_gamg->Nlevels; level++) {
8279566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Parr[level]));
8289566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Aarr[level]));
8295b89ad90SMark F. Adams     }
830806fa848SBarry Smith   } else {
8315f8cf99dSMark F. Adams     KSP smoother;
8320ed2132dSStefano Zampini 
83363a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix));
8349566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
8359566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
8369566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPPREONLY));
8379566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8385f8cf99dSMark F. Adams   }
839849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8415b89ad90SMark F. Adams }
8425b89ad90SMark F. Adams 
8435b89ad90SMark F. Adams /*
8445b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8455b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8465b89ad90SMark F. Adams 
8475b89ad90SMark F. Adams    Input Parameter:
8485b89ad90SMark F. Adams .  pc - the preconditioner context
8495b89ad90SMark F. Adams 
8505b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8515b89ad90SMark F. Adams */
852d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_GAMG(PC pc)
853d71ae5a4SJacob Faibussowitsch {
8545b89ad90SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
8555b89ad90SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
8565b89ad90SMark F. Adams 
8575b89ad90SMark F. Adams   PetscFunctionBegin;
8589566063dSJacob Faibussowitsch   PetscCall(PCReset_GAMG(pc));
8591baa6e33SBarry Smith   if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc));
8609566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->ops));
8619566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
8629566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg));
8632e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL));
8642e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL));
8652e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL));
8662e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL));
8672e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL));
868e1cf1444SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", NULL));
8692e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL));
8702e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL));
871d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", NULL));
8722e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL));
8732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL));
8742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL));
8752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL));
8762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL));
8772e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL));
8782e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL));
8792e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL));
8809566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
8813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8825b89ad90SMark F. Adams }
8835b89ad90SMark F. Adams 
884676e1743SMark F. Adams /*@
885f1580f4eSBarry Smith   PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG`
886676e1743SMark F. Adams 
887c3339decSBarry Smith   Logically Collective
888676e1743SMark F. Adams 
889676e1743SMark F. Adams   Input Parameters:
8901cc46a46SBarry Smith + pc - the preconditioner context
8911cc46a46SBarry Smith - n  - the number of equations
892676e1743SMark F. Adams 
893676e1743SMark F. Adams   Options Database Key:
894147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit
895676e1743SMark F. Adams 
89620f4b53cSBarry Smith   Level: intermediate
89720f4b53cSBarry Smith 
898f1580f4eSBarry Smith   Note:
899f1580f4eSBarry Smith   `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
900cab9ed1eSBarry Smith   that has degrees of freedom
901cab9ed1eSBarry Smith 
902*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
903676e1743SMark F. Adams @*/
904d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
905d71ae5a4SJacob Faibussowitsch {
906676e1743SMark F. Adams   PetscFunctionBegin;
907676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
908cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n));
9093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
910676e1743SMark F. Adams }
911676e1743SMark F. Adams 
912d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
913d71ae5a4SJacob Faibussowitsch {
914c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
915c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
916676e1743SMark F. Adams 
917676e1743SMark F. Adams   PetscFunctionBegin;
9189d5b6da9SMark F. Adams   if (n > 0) pc_gamg->min_eq_proc = n;
9193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
920676e1743SMark F. Adams }
921676e1743SMark F. Adams 
922389730f3SMark F. Adams /*@
923f1580f4eSBarry Smith   PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG`
924389730f3SMark F. Adams 
925c3339decSBarry Smith   Collective
926389730f3SMark F. Adams 
927389730f3SMark F. Adams   Input Parameters:
9281cc46a46SBarry Smith + pc - the preconditioner context
9291cc46a46SBarry Smith - n  - maximum number of equations to aim for
930389730f3SMark F. Adams 
931389730f3SMark F. Adams   Options Database Key:
932147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit
933389730f3SMark F. Adams 
93420f4b53cSBarry Smith   Level: intermediate
93520f4b53cSBarry Smith 
936f1580f4eSBarry Smith   Note:
937f1580f4eSBarry Smith   For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
93874329af1SBarry Smith   has less than 1000 unknowns.
93974329af1SBarry Smith 
940*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
941389730f3SMark F. Adams @*/
942d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
943d71ae5a4SJacob Faibussowitsch {
944389730f3SMark F. Adams   PetscFunctionBegin;
945389730f3SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
946cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n));
9473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
948389730f3SMark F. Adams }
949389730f3SMark F. Adams 
950d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
951d71ae5a4SJacob Faibussowitsch {
952389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
953389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
954389730f3SMark F. Adams 
955389730f3SMark F. Adams   PetscFunctionBegin;
9569d5b6da9SMark F. Adams   if (n > 0) pc_gamg->coarse_eq_limit = n;
9573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
958389730f3SMark F. Adams }
959389730f3SMark F. Adams 
960676e1743SMark F. Adams /*@
961f1580f4eSBarry Smith   PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use
962676e1743SMark F. Adams 
963c3339decSBarry Smith   Collective
964676e1743SMark F. Adams 
965676e1743SMark F. Adams   Input Parameters:
9661cc46a46SBarry Smith + pc - the preconditioner context
967f1580f4eSBarry Smith - n  - `PETSC_TRUE` or `PETSC_FALSE`
968676e1743SMark F. Adams 
969676e1743SMark F. Adams   Options Database Key:
970147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning
971676e1743SMark F. Adams 
97220f4b53cSBarry Smith   Level: intermediate
97320f4b53cSBarry Smith 
974f1580f4eSBarry Smith   Note:
975f1580f4eSBarry Smith   This will generally improve the loading balancing of the work on each level
976cab9ed1eSBarry Smith 
977*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
978676e1743SMark F. Adams @*/
979d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
980d71ae5a4SJacob Faibussowitsch {
981676e1743SMark F. Adams   PetscFunctionBegin;
982676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
983cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n));
9843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
985676e1743SMark F. Adams }
986676e1743SMark F. Adams 
987d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
988d71ae5a4SJacob Faibussowitsch {
989c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
990c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
991676e1743SMark F. Adams 
992676e1743SMark F. Adams   PetscFunctionBegin;
9939d5b6da9SMark F. Adams   pc_gamg->repart = n;
9943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
995676e1743SMark F. Adams }
996676e1743SMark F. Adams 
997dfd5c07aSMark F. Adams /*@
998f1580f4eSBarry Smith   PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process
99918c3aa7eSMark 
1000c3339decSBarry Smith   Collective
100118c3aa7eSMark 
100218c3aa7eSMark   Input Parameters:
100318c3aa7eSMark + pc - the preconditioner context
1004e1cf1444SMark Adams - b  - flag
100518c3aa7eSMark 
100618c3aa7eSMark   Options Database Key:
1007147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
100818c3aa7eSMark 
100920f4b53cSBarry Smith   Level: advanced
101020f4b53cSBarry Smith 
101118c3aa7eSMark   Notes:
10127e6512fdSJed 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$.
1013f1580f4eSBarry Smith   Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
1014f1580f4eSBarry Smith   If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process
1015f1580f4eSBarry Smith   This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used.
1016efe053fcSJed Brown   It became default in PETSc 3.17.
101718c3aa7eSMark 
1018*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `PCGAMGSetRecomputeEstEig()`
101918c3aa7eSMark @*/
1020e1cf1444SMark Adams PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool b)
1021d71ae5a4SJacob Faibussowitsch {
102218c3aa7eSMark   PetscFunctionBegin;
102318c3aa7eSMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1024e1cf1444SMark Adams   PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, b));
10253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102618c3aa7eSMark }
102718c3aa7eSMark 
1028e1cf1444SMark Adams static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool b)
1029d71ae5a4SJacob Faibussowitsch {
103018c3aa7eSMark   PC_MG   *mg      = (PC_MG *)pc->data;
103118c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
103218c3aa7eSMark 
103318c3aa7eSMark   PetscFunctionBegin;
1034e1cf1444SMark Adams   pc_gamg->use_sa_esteig = b;
1035e1cf1444SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1036e1cf1444SMark Adams }
1037e1cf1444SMark Adams 
1038e1cf1444SMark Adams /*@
1039e1cf1444SMark Adams   PCGAMGSetRecomputeEstEig - Set flag for Chebyshev smoothers to recompute the eigen estimates
1040e1cf1444SMark Adams 
1041e1cf1444SMark Adams   Collective
1042e1cf1444SMark Adams 
1043e1cf1444SMark Adams   Input Parameters:
1044e1cf1444SMark Adams + pc - the preconditioner context
1045e1cf1444SMark Adams - b  - flag
1046e1cf1444SMark Adams 
1047e1cf1444SMark Adams   Options Database Key:
1048e1cf1444SMark Adams . -pc_gamg_recompute_esteig <true> - use the eigen estimate
1049e1cf1444SMark Adams 
1050e1cf1444SMark Adams   Level: advanced
1051e1cf1444SMark Adams 
1052*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
1053e1cf1444SMark Adams @*/
1054e1cf1444SMark Adams PetscErrorCode PCGAMGSetRecomputeEstEig(PC pc, PetscBool b)
1055e1cf1444SMark Adams {
1056e1cf1444SMark Adams   PetscFunctionBegin;
1057e1cf1444SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1058e1cf1444SMark Adams   PetscTryMethod(pc, "PCGAMGSetRecomputeEstEig_C", (PC, PetscBool), (pc, b));
1059e1cf1444SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1060e1cf1444SMark Adams }
1061e1cf1444SMark Adams 
1062e1cf1444SMark Adams static PetscErrorCode PCGAMGSetRecomputeEstEig_GAMG(PC pc, PetscBool b)
1063e1cf1444SMark Adams {
1064e1cf1444SMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1065e1cf1444SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1066e1cf1444SMark Adams 
1067e1cf1444SMark Adams   PetscFunctionBegin;
1068e1cf1444SMark Adams   pc_gamg->recompute_esteig = b;
10693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107018c3aa7eSMark }
107118c3aa7eSMark 
107218c3aa7eSMark /*@
1073f1580f4eSBarry Smith   PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY?
107418c3aa7eSMark 
1075c3339decSBarry Smith   Collective
107618c3aa7eSMark 
107718c3aa7eSMark   Input Parameters:
107818c3aa7eSMark + pc   - the preconditioner context
1079feefa0e1SJacob Faibussowitsch . emax - max eigenvalue
1080feefa0e1SJacob Faibussowitsch - emin - min eigenvalue
108118c3aa7eSMark 
108218c3aa7eSMark   Options Database Key:
1083147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
108418c3aa7eSMark 
108518c3aa7eSMark   Level: intermediate
108618c3aa7eSMark 
1087*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetUseSAEstEig()`
108818c3aa7eSMark @*/
1089d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin)
1090d71ae5a4SJacob Faibussowitsch {
109118c3aa7eSMark   PetscFunctionBegin;
109218c3aa7eSMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1093cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin));
10943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109518c3aa7eSMark }
109641ffd417SStefano Zampini 
1097d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin)
1098d71ae5a4SJacob Faibussowitsch {
109918c3aa7eSMark   PC_MG   *mg      = (PC_MG *)pc->data;
110018c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
110118c3aa7eSMark 
110218c3aa7eSMark   PetscFunctionBegin;
11035f80ce2aSJacob 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);
11045f80ce2aSJacob 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);
110518c3aa7eSMark   pc_gamg->emax = emax;
110618c3aa7eSMark   pc_gamg->emin = emin;
11073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110818c3aa7eSMark }
110918c3aa7eSMark 
111018c3aa7eSMark /*@
1111f1580f4eSBarry Smith   PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner
1112dfd5c07aSMark F. Adams 
1113c3339decSBarry Smith   Collective
1114dfd5c07aSMark F. Adams 
1115dfd5c07aSMark F. Adams   Input Parameters:
11161cc46a46SBarry Smith + pc - the preconditioner context
1117f1580f4eSBarry Smith - n  - `PETSC_TRUE` or `PETSC_FALSE`
1118dfd5c07aSMark F. Adams 
1119dfd5c07aSMark F. Adams   Options Database Key:
1120147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1121dfd5c07aSMark F. Adams 
1122dfd5c07aSMark F. Adams   Level: intermediate
1123dfd5c07aSMark F. Adams 
1124f1580f4eSBarry Smith   Note:
1125147403d9SBarry Smith   May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1126cab9ed1eSBarry Smith   rebuilding the preconditioner quicker.
1127cab9ed1eSBarry Smith 
1128*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`
1129dfd5c07aSMark F. Adams @*/
1130d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1131d71ae5a4SJacob Faibussowitsch {
1132dfd5c07aSMark F. Adams   PetscFunctionBegin;
1133dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1134cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n));
11353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1136dfd5c07aSMark F. Adams }
1137dfd5c07aSMark F. Adams 
1138d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1139d71ae5a4SJacob Faibussowitsch {
1140dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1141dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1142dfd5c07aSMark F. Adams 
1143dfd5c07aSMark F. Adams   PetscFunctionBegin;
1144dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
11453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1146dfd5c07aSMark F. Adams }
1147dfd5c07aSMark F. Adams 
1148ffc955d6SMark F. Adams /*@
1149f1580f4eSBarry 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
1150f1580f4eSBarry Smith   used as the smoother
1151ffc955d6SMark F. Adams 
1152c3339decSBarry Smith   Collective
1153ffc955d6SMark F. Adams 
1154ffc955d6SMark F. Adams   Input Parameters:
1155cab9ed1eSBarry Smith + pc  - the preconditioner context
1156f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not
1157ffc955d6SMark F. Adams 
1158ffc955d6SMark F. Adams   Options Database Key:
1159147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1160ffc955d6SMark F. Adams 
1161ffc955d6SMark F. Adams   Level: intermediate
1162ffc955d6SMark F. Adams 
1163*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCASM`, `PCSetType`
1164ffc955d6SMark F. Adams @*/
1165d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1166d71ae5a4SJacob Faibussowitsch {
1167ffc955d6SMark F. Adams   PetscFunctionBegin;
1168ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1169cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg));
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1171ffc955d6SMark F. Adams }
1172ffc955d6SMark F. Adams 
1173d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1174d71ae5a4SJacob Faibussowitsch {
1175ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1176ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1177ffc955d6SMark F. Adams 
1178ffc955d6SMark F. Adams   PetscFunctionBegin;
1179cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
11803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1181ffc955d6SMark F. Adams }
1182ffc955d6SMark F. Adams 
1183171cca9aSMark Adams /*@
1184d529f056Smarkadams4   PCGAMGSetParallelCoarseGridSolve - allow a parallel coarse grid solver
1185171cca9aSMark Adams 
1186c3339decSBarry Smith   Collective
1187171cca9aSMark Adams 
1188171cca9aSMark Adams   Input Parameters:
1189171cca9aSMark Adams + pc  - the preconditioner context
1190f1580f4eSBarry Smith - flg - `PETSC_TRUE` to not force coarse grid onto one processor
1191171cca9aSMark Adams 
1192171cca9aSMark Adams   Options Database Key:
1193d529f056Smarkadams4 . -pc_gamg_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1194171cca9aSMark Adams 
1195171cca9aSMark Adams   Level: intermediate
1196171cca9aSMark Adams 
1197*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`
1198171cca9aSMark Adams @*/
1199d529f056Smarkadams4 PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC pc, PetscBool flg)
1200d71ae5a4SJacob Faibussowitsch {
1201171cca9aSMark Adams   PetscFunctionBegin;
1202171cca9aSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1203d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGSetParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg));
12043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1205171cca9aSMark Adams }
1206171cca9aSMark Adams 
1207d529f056Smarkadams4 static PetscErrorCode PCGAMGSetParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1208d71ae5a4SJacob Faibussowitsch {
1209171cca9aSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1210171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1211171cca9aSMark Adams 
1212171cca9aSMark Adams   PetscFunctionBegin;
1213171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
12143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1215ffc955d6SMark F. Adams }
1216ffc955d6SMark F. Adams 
12174ef23d27SMark F. Adams /*@
1218f1580f4eSBarry Smith   PCGAMGSetCpuPinCoarseGrids - pin the coarse grids created in `PCGAMG` to run only on the CPU since the problems may be too small to run efficiently on the GPUs
1219ce7c7f2fSMark Adams 
1220c3339decSBarry Smith   Collective
1221ce7c7f2fSMark Adams 
1222ce7c7f2fSMark Adams   Input Parameters:
1223ce7c7f2fSMark Adams + pc  - the preconditioner context
1224f1580f4eSBarry Smith - flg - `PETSC_TRUE` to pin coarse grids to the CPU
1225ce7c7f2fSMark Adams 
1226ce7c7f2fSMark Adams   Options Database Key:
1227147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1228ce7c7f2fSMark Adams 
1229ce7c7f2fSMark Adams   Level: intermediate
1230ce7c7f2fSMark Adams 
1231*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetParallelCoarseGridSolve()`
1232ce7c7f2fSMark Adams @*/
1233d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1234d71ae5a4SJacob Faibussowitsch {
1235ce7c7f2fSMark Adams   PetscFunctionBegin;
1236ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1237cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg));
12383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1239ce7c7f2fSMark Adams }
1240ce7c7f2fSMark Adams 
1241d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1242d71ae5a4SJacob Faibussowitsch {
1243ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1244ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1245ce7c7f2fSMark Adams 
1246ce7c7f2fSMark Adams   PetscFunctionBegin;
1247ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
12483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1249ce7c7f2fSMark Adams }
1250ce7c7f2fSMark Adams 
1251ce7c7f2fSMark Adams /*@
1252147403d9SBarry Smith   PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1253ce7c7f2fSMark Adams 
1254c3339decSBarry Smith   Collective
1255ce7c7f2fSMark Adams 
1256ce7c7f2fSMark Adams   Input Parameters:
1257ce7c7f2fSMark Adams + pc  - the preconditioner context
1258f1580f4eSBarry Smith - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD`
1259ce7c7f2fSMark Adams 
1260ce7c7f2fSMark Adams   Options Database Key:
1261147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1262ce7c7f2fSMark Adams 
1263ce7c7f2fSMark Adams   Level: intermediate
1264ce7c7f2fSMark Adams 
1265*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD`
1266ce7c7f2fSMark Adams @*/
1267d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1268d71ae5a4SJacob Faibussowitsch {
1269ce7c7f2fSMark Adams   PetscFunctionBegin;
1270ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1271cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg));
12723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1273ce7c7f2fSMark Adams }
1274ce7c7f2fSMark Adams 
1275d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1276d71ae5a4SJacob Faibussowitsch {
1277ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1278ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1279ce7c7f2fSMark Adams 
1280ce7c7f2fSMark Adams   PetscFunctionBegin;
1281ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
12823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1283ce7c7f2fSMark Adams }
1284ce7c7f2fSMark Adams 
1285ce7c7f2fSMark Adams /*@
1286f1580f4eSBarry Smith   PCGAMGSetNlevels -  Sets the maximum number of levels `PCGAMG` will use
12874ef23d27SMark F. Adams 
128820f4b53cSBarry Smith   Not Collective
12894ef23d27SMark F. Adams 
12904ef23d27SMark F. Adams   Input Parameters:
12911cc46a46SBarry Smith + pc - the preconditioner
12921cc46a46SBarry Smith - n  - the maximum number of levels to use
12934ef23d27SMark F. Adams 
12944ef23d27SMark F. Adams   Options Database Key:
1295147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow
12964ef23d27SMark F. Adams 
12974ef23d27SMark F. Adams   Level: intermediate
12984ef23d27SMark F. Adams 
1299feefa0e1SJacob Faibussowitsch   Developer Notes:
1300f1580f4eSBarry Smith   Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`
1301f1580f4eSBarry Smith 
1302*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`
13034ef23d27SMark F. Adams @*/
1304d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1305d71ae5a4SJacob Faibussowitsch {
13064ef23d27SMark F. Adams   PetscFunctionBegin;
13074ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1308cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n));
13093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13104ef23d27SMark F. Adams }
13114ef23d27SMark F. Adams 
1312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1313d71ae5a4SJacob Faibussowitsch {
13144ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
13154ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
13164ef23d27SMark F. Adams 
13174ef23d27SMark F. Adams   PetscFunctionBegin;
13189d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
13193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13204ef23d27SMark F. Adams }
13214ef23d27SMark F. Adams 
13223542efc5SMark F. Adams /*@
13233542efc5SMark F. Adams   PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
13243542efc5SMark F. Adams 
132520f4b53cSBarry Smith   Not Collective
13263542efc5SMark F. Adams 
13273542efc5SMark F. Adams   Input Parameters:
13281cc46a46SBarry Smith + pc - the preconditioner context
1329c9567895SMark . 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
1330055c8bd0SJed Brown - n  - number of threshold values provided in array
13313542efc5SMark F. Adams 
13323542efc5SMark F. Adams   Options Database Key:
1333147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges
13343542efc5SMark F. Adams 
133520f4b53cSBarry Smith   Level: intermediate
133620f4b53cSBarry Smith 
133795452b02SPatrick Sanan   Notes:
1338af3c827dSMark 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.
1339f1580f4eSBarry Smith   Before coarsening or aggregating the graph, `PCGAMG` 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.
1340cab9ed1eSBarry Smith 
134120f4b53cSBarry Smith   If `n` is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening.
1342f1580f4eSBarry Smith   In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`.
134320f4b53cSBarry Smith   If `n` is greater than the total number of levels, the excess entries in threshold will not be used.
13443542efc5SMark F. Adams 
1345*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetMinDegreeOrderingMISk()`, `PCGAMGSetThresholdScale()`
13463542efc5SMark F. Adams @*/
1347d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1348d71ae5a4SJacob Faibussowitsch {
13493542efc5SMark F. Adams   PetscFunctionBegin;
13503542efc5SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13514f572ea9SToby Isaac   if (n) PetscAssertPointer(v, 2);
1352cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n));
13533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13543542efc5SMark F. Adams }
13553542efc5SMark F. Adams 
1356d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1357d71ae5a4SJacob Faibussowitsch {
1358c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1359c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1360c1eae691SMark Adams   PetscInt i;
1361c1eae691SMark Adams   PetscFunctionBegin;
1362055c8bd0SJed Brown   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1363055c8bd0SJed Brown   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
13643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1365c1eae691SMark Adams }
1366c1eae691SMark Adams 
1367c1eae691SMark Adams /*@
1368f1580f4eSBarry Smith   PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids
1369c9567895SMark 
1370c3339decSBarry Smith   Collective
1371c9567895SMark 
1372c9567895SMark   Input Parameters:
1373c9567895SMark + pc - the preconditioner context
1374f1580f4eSBarry Smith . v  - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA
1375c9567895SMark - n  - number of values provided in array
1376c9567895SMark 
1377c9567895SMark   Options Database Key:
1378147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1379c9567895SMark 
1380c9567895SMark   Level: intermediate
1381c9567895SMark 
1382*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`
1383c9567895SMark @*/
1384d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1385d71ae5a4SJacob Faibussowitsch {
1386c9567895SMark   PetscFunctionBegin;
1387c9567895SMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13884f572ea9SToby Isaac   if (n) PetscAssertPointer(v, 2);
1389cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n));
13903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1391c9567895SMark }
1392c9567895SMark 
1393d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1394d71ae5a4SJacob Faibussowitsch {
1395c9567895SMark   PC_MG   *mg      = (PC_MG *)pc->data;
1396c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1397c9567895SMark   PetscInt i;
1398c9567895SMark   PetscFunctionBegin;
1399c9567895SMark   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1400c9567895SMark   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
14013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1402c9567895SMark }
1403c9567895SMark 
1404c9567895SMark /*@
1405c1eae691SMark Adams   PCGAMGSetThresholdScale - Relative threshold reduction at each level
1406c1eae691SMark Adams 
140720f4b53cSBarry Smith   Not Collective
1408c1eae691SMark Adams 
1409c1eae691SMark Adams   Input Parameters:
1410c1eae691SMark Adams + pc - the preconditioner context
1411feefa0e1SJacob Faibussowitsch - v  - the threshold value reduction, usually < 1.0
1412c1eae691SMark Adams 
1413c1eae691SMark Adams   Options Database Key:
1414147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1415c1eae691SMark Adams 
141620f4b53cSBarry Smith   Level: advanced
141720f4b53cSBarry Smith 
1418f1580f4eSBarry Smith   Note:
1419f1580f4eSBarry Smith   The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`.
1420f1580f4eSBarry Smith   This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`.
1421055c8bd0SJed Brown 
1422*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`
1423c1eae691SMark Adams @*/
1424d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1425d71ae5a4SJacob Faibussowitsch {
14263542efc5SMark F. Adams   PetscFunctionBegin;
1427c1eae691SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1428cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v));
14293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1430c1eae691SMark Adams }
1431c1eae691SMark Adams 
1432d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1433d71ae5a4SJacob Faibussowitsch {
1434c1eae691SMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1435c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1436c1eae691SMark Adams   PetscFunctionBegin;
1437c1eae691SMark Adams   pc_gamg->threshold_scale = v;
14383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14393542efc5SMark F. Adams }
14403542efc5SMark F. Adams 
1441e20c40e8SBarry Smith /*@C
1442f1580f4eSBarry Smith   PCGAMGSetType - Set the type of algorithm `PCGAMG` should use
1443676e1743SMark F. Adams 
1444c3339decSBarry Smith   Collective
1445676e1743SMark F. Adams 
1446676e1743SMark F. Adams   Input Parameters:
1447c60c7ad4SBarry Smith + pc   - the preconditioner context
1448f1580f4eSBarry Smith - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL`
1449676e1743SMark F. Adams 
1450676e1743SMark F. Adams   Options Database Key:
14515e4ac8c8Smarkadams4 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply - only agg supported
1452676e1743SMark F. Adams 
1453676e1743SMark F. Adams   Level: intermediate
1454676e1743SMark F. Adams 
1455*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType`
1456676e1743SMark F. Adams @*/
1457d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1458d71ae5a4SJacob Faibussowitsch {
1459676e1743SMark F. Adams   PetscFunctionBegin;
1460676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1461cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type));
14623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1463676e1743SMark F. Adams }
1464676e1743SMark F. Adams 
1465e20c40e8SBarry Smith /*@C
1466f1580f4eSBarry Smith   PCGAMGGetType - Get the type of algorithm `PCGAMG` will use
1467c60c7ad4SBarry Smith 
1468c3339decSBarry Smith   Collective
1469c60c7ad4SBarry Smith 
1470c60c7ad4SBarry Smith   Input Parameter:
1471c60c7ad4SBarry Smith . pc - the preconditioner context
1472c60c7ad4SBarry Smith 
1473c60c7ad4SBarry Smith   Output Parameter:
1474c60c7ad4SBarry Smith . type - the type of algorithm used
1475c60c7ad4SBarry Smith 
1476c60c7ad4SBarry Smith   Level: intermediate
1477c60c7ad4SBarry Smith 
1478*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType`
1479c60c7ad4SBarry Smith @*/
1480d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1481d71ae5a4SJacob Faibussowitsch {
1482c60c7ad4SBarry Smith   PetscFunctionBegin;
1483c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1484cac4c232SBarry Smith   PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type));
14853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1486c60c7ad4SBarry Smith }
1487c60c7ad4SBarry Smith 
1488d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1489d71ae5a4SJacob Faibussowitsch {
1490c60c7ad4SBarry Smith   PC_MG   *mg      = (PC_MG *)pc->data;
1491c60c7ad4SBarry Smith   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1492c60c7ad4SBarry Smith 
1493c60c7ad4SBarry Smith   PetscFunctionBegin;
1494c60c7ad4SBarry Smith   *type = pc_gamg->type;
14953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1496c60c7ad4SBarry Smith }
1497c60c7ad4SBarry Smith 
1498d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1499d71ae5a4SJacob Faibussowitsch {
15001ab5ffc9SJed Brown   PC_MG   *mg      = (PC_MG *)pc->data;
15011ab5ffc9SJed Brown   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
15025f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PC);
1503676e1743SMark F. Adams 
1504676e1743SMark F. Adams   PetscFunctionBegin;
1505c60c7ad4SBarry Smith   pc_gamg->type = type;
15069566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(GAMGList, type, &r));
15076adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type);
15081ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
15099566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
15109566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps)));
1511e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1512da81f932SPierre Jolivet     /* cleaning up common data in pc_gamg - this should disappear someday */
15133ae0bb68SMark Adams     pc_gamg->data_cell_cols      = 0;
15143ae0bb68SMark Adams     pc_gamg->data_cell_rows      = 0;
15153ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
15163ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
15179566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
15183ae0bb68SMark Adams     pc_gamg->data_sz = 0;
15191ab5ffc9SJed Brown   }
15209566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
15219566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name));
15229566063dSJacob Faibussowitsch   PetscCall((*r)(pc));
15233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1524676e1743SMark F. Adams }
1525676e1743SMark F. Adams 
1526d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer)
1527d71ae5a4SJacob Faibussowitsch {
15285adeb434SBarry Smith   PC_MG    *mg      = (PC_MG *)pc->data;
15295adeb434SBarry Smith   PC_GAMG  *pc_gamg = (PC_GAMG *)mg->innerctx;
1530e7d4b4cbSMark Adams   PetscReal gc = 0, oc = 0;
153190db8557SMark Adams 
15325adeb434SBarry Smith   PetscFunctionBegin;
15339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "    GAMG specific options\n"));
15349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold for dropping small values in graph on each level ="));
15359566063dSJacob Faibussowitsch   for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i]));
15369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
15379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale));
153848a46eb9SPierre Jolivet   if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates from coarsening process to define subdomains for PCASM\n"));
153948a46eb9SPierre Jolivet   if (pc_gamg->use_parallel_coarse_grid_solver) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n"));
15401baa6e33SBarry Smith   if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer));
15419566063dSJacob Faibussowitsch   PetscCall(PCMGGetGridComplexity(pc, &gc, &oc));
154263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Complexity:    grid = %g    operator = %g\n", (double)gc, (double)oc));
15433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15445adeb434SBarry Smith }
15455adeb434SBarry Smith 
154666976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject)
1547d71ae5a4SJacob Faibussowitsch {
1548676e1743SMark F. Adams   PC_MG             *mg      = (PC_MG *)pc->data;
1549676e1743SMark F. Adams   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
15507e6512fdSJed Brown   PetscBool          flag;
15513b4367a7SBarry Smith   MPI_Comm           comm;
155218c3aa7eSMark   char               prefix[256], tname[32];
1553c1eae691SMark Adams   PetscInt           i, n;
155414a9496bSBarry Smith   const char        *pcpre;
15550a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL};
15565b89ad90SMark F. Adams   PetscFunctionBegin;
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1558d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options");
15595e4ac8c8Smarkadams4   PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method (only 'agg' supported and useful)", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
15601baa6e33SBarry Smith   if (flag) PetscCall(PCGAMGSetType(pc, tname));
15619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL));
1562f1580f4eSBarry Smith   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));
1563e1cf1444SMark Adams   PetscCall(PetscOptionsBool("-pc_gamg_recompute_esteig", "Set flag to recompute eigen estimates for Chebyshev when matrix changes", "PCGAMGSetRecomputeEstEig", pc_gamg->recompute_esteig, &pc_gamg->recompute_esteig, NULL));
15649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL));
15659566063dSJacob 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));
1566d529f056Smarkadams4   PetscCall(PetscOptionsBool("-pc_gamg_parallel_coarse_grid_solver", "Use parallel coarse grid solver (otherwise put last grid on one process)", "PCGAMGSetParallelCoarseGridSolve", pc_gamg->use_parallel_coarse_grid_solver, &pc_gamg->use_parallel_coarse_grid_solver, NULL));
15679566063dSJacob 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));
15689371c9d4SSatish Balay   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,
15699371c9d4SSatish Balay                              (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL));
15709566063dSJacob 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));
15719566063dSJacob 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));
15729566063dSJacob 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));
157318c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
15749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag));
157518c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1576efd3c5ceSMark Adams     if (!flag) n = 1;
1577c1eae691SMark Adams     i = n;
1578d71ae5a4SJacob Faibussowitsch     do {
1579d71ae5a4SJacob Faibussowitsch       pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1580d71ae5a4SJacob Faibussowitsch     } while (++i < PETSC_MG_MAXLEVELS);
1581c1eae691SMark Adams   }
1582c9567895SMark   n = PETSC_MG_MAXLEVELS;
15839566063dSJacob 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));
1584c9567895SMark   if (!flag) i = 0;
1585c9567895SMark   else i = n;
1586d71ae5a4SJacob Faibussowitsch   do {
1587d71ae5a4SJacob Faibussowitsch     pc_gamg->level_reduction_factors[i] = -1;
1588d71ae5a4SJacob Faibussowitsch   } while (++i < PETSC_MG_MAXLEVELS);
15899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL));
159018c3aa7eSMark   {
159118c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
159218c3aa7eSMark     n                    = 2;
15939566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag));
159418c3aa7eSMark     if (flag) {
159508401ef6SPierre Jolivet       PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
15969566063dSJacob Faibussowitsch       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
159718c3aa7eSMark     }
159818c3aa7eSMark   }
1599b7cbab4eSMark Adams   /* set options for subtype */
1600dbbe0bcdSBarry Smith   PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject));
160118c3aa7eSMark 
16029566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
16039566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
1604d0609cedSBarry Smith   PetscOptionsHeadEnd();
16053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16065b89ad90SMark F. Adams }
16075b89ad90SMark F. Adams 
16085b89ad90SMark F. Adams /*MC
16091cc46a46SBarry Smith   PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
16105b89ad90SMark F. Adams 
1611280d9858SJed Brown   Options Database Keys:
16125e4ac8c8Smarkadams4 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical (only smoothed aggregation, agg, supported)
1613da81f932SPierre Jolivet . -pc_gamg_repartition  <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined
161421d928e4Smarkadams4 . -pc_gamg_asm_use_agg <bool,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1615da81f932SPierre Jolivet . -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI ranks used directly on the coarse grids so that there are around <limit>
1616cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
16172d776b49SBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
161821d928e4Smarkadams4 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true)
161921d928e4Smarkadams4 . -pc_gamg_threshold[] <thresh,default=[-1,...]> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 does no filtering)
16202d776b49SBarry Smith - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1621cab9ed1eSBarry Smith 
1622f1580f4eSBarry Smith   Options Database Keys for Aggregation:
1623cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1624d529f056Smarkadams4 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1625d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1626d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1627d529f056Smarkadams4 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1628cab9ed1eSBarry Smith 
1629f1580f4eSBarry Smith   Options Database Keys for Multigrid:
1630a9f5add0SYANG Zongze + -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()`
1631db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1632db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
163321d928e4Smarkadams4 - -pc_mg_levels <levels> - Number of levels of multigrid to use. GAMG has a heuristic so pc_mg_levels is not usually used with GAMG
16345b89ad90SMark F. Adams 
163520f4b53cSBarry Smith   Level: intermediate
163620f4b53cSBarry Smith 
163795452b02SPatrick Sanan   Notes:
1638f1580f4eSBarry Smith   To obtain good performance for `PCGAMG` for vector valued problems you must
1639f1580f4eSBarry Smith   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point
1640f1580f4eSBarry Smith   call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1641f1580f4eSBarry Smith 
164204c3f3b8SBarry Smith   The many options for `PCMG` also work directly for `PCGAMG` such as controlling the smoothers on each level etc.
164304c3f3b8SBarry Smith 
164404c3f3b8SBarry Smith   See [the Users Manual section on PCGAMG](sec_amg) and [the Users Manual section on PCMG](sec_mg)for more details.
16451cc46a46SBarry Smith 
1646*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`,
1647d529f056Smarkadams4           `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
16485b89ad90SMark F. Adams M*/
1649d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1650d71ae5a4SJacob Faibussowitsch {
16515b89ad90SMark F. Adams   PC_GAMG *pc_gamg;
16525b89ad90SMark F. Adams   PC_MG   *mg;
16535b89ad90SMark F. Adams 
16545b89ad90SMark F. Adams   PetscFunctionBegin;
16551c1aac46SBarry Smith   /* register AMG type */
16569566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
16571c1aac46SBarry Smith 
16585b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
16599566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
16609566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
16615b89ad90SMark F. Adams 
16625b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
16634dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg));
16649566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
16655b89ad90SMark F. Adams   mg           = (PC_MG *)pc->data;
16665b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
16675b89ad90SMark F. Adams 
16684dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg->ops));
16691ab5ffc9SJed Brown 
16709d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
16719d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
16720a545947SLisandro Dalcin   pc_gamg->data    = NULL;
16735b89ad90SMark F. Adams 
16749d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
16755b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
16765b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
16775b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
16785b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
16795adeb434SBarry Smith   mg->view                = PCView_GAMG;
16805b89ad90SMark F. Adams 
16819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG));
16829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG));
16839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG));
16849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG));
16859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG));
1686e1cf1444SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", PCGAMGSetRecomputeEstEig_GAMG));
16879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG));
16889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG));
1689d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", PCGAMGSetParallelCoarseGridSolve_GAMG));
16909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG));
16919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG));
16929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG));
16939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG));
16949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG));
16959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG));
16969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG));
16979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG));
16989d5b6da9SMark F. Adams   pc_gamg->repart                          = PETSC_FALSE;
169921d928e4Smarkadams4   pc_gamg->reuse_prol                      = PETSC_TRUE;
17000c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm                 = PETSC_FALSE;
1701171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1702a0095786SMark   pc_gamg->cpu_pin_coarse_grids            = PETSC_FALSE;
1703a0095786SMark   pc_gamg->layout_type                     = PCGAMG_LAYOUT_SPREAD;
1704038f3aa4SMark F. Adams   pc_gamg->min_eq_proc                     = 50;
170525a145a7SMark Adams   pc_gamg->coarse_eq_limit                 = 50;
170653134ebeSMark Adams   for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1;
1707c1eae691SMark Adams   pc_gamg->threshold_scale  = 1.;
170818c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
17099ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
17107e6512fdSJed Brown   pc_gamg->use_sa_esteig    = PETSC_TRUE;
1711e1cf1444SMark Adams   pc_gamg->recompute_esteig = PETSC_TRUE;
171218c3aa7eSMark   pc_gamg->emin             = 0;
171318c3aa7eSMark   pc_gamg->emax             = 0;
171418c3aa7eSMark 
1715c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
17169d5b6da9SMark F. Adams 
1717bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
17189566063dSJacob Faibussowitsch   PetscCall(PCGAMGSetType(pc, PCGAMGAGG));
17193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17205b89ad90SMark F. Adams }
17213e3471ccSMark Adams 
17223e3471ccSMark Adams /*@C
1723f1580f4eSBarry Smith   PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called
1724f1580f4eSBarry Smith   from `PCInitializePackage()`.
17253e3471ccSMark Adams 
17263e3471ccSMark Adams   Level: developer
17273e3471ccSMark Adams 
1728*562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscInitialize()`
17293e3471ccSMark Adams @*/
1730d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGInitializePackage(void)
1731d71ae5a4SJacob Faibussowitsch {
17324555aa8cSStefano Zampini   PetscInt l;
17333e3471ccSMark Adams 
17343e3471ccSMark Adams   PetscFunctionBegin;
17353ba16761SJacob Faibussowitsch   if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
17363e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
17379566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO));
17389566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG));
17399566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical));
17409566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1741c1c463dbSMark Adams 
1742c1c463dbSMark Adams   /* general events */
1743849bee69SMark Adams   PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]));
1744849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]));
1745849bee69SMark Adams   PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]));
1746849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]));
1747849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]));
1748849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]));
1749849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]));
1750849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]));
1751849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]));
1752849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]));
1753849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]));
1754849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]));
1755849bee69SMark Adams   PetscCall(PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]));
1756849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */
1757849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */
1758849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */
17594555aa8cSStefano Zampini   for (l = 0; l < PETSC_MG_MAXLEVELS; l++) {
17604555aa8cSStefano Zampini     char ename[32];
17615b89ad90SMark F. Adams 
176263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l));
17639566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
176463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l));
17659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
176663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l));
17679566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
17684555aa8cSStefano Zampini   }
17694555aa8cSStefano Zampini #if defined(GAMG_STAGES)
1770849bee69SMark Adams   { /* create timer stages */
17715b89ad90SMark F. Adams     char str[32];
1772a364092eSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0));
17739566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
17745b89ad90SMark F. Adams   }
17755b89ad90SMark F. Adams #endif
17763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17773e3471ccSMark Adams }
17783e3471ccSMark Adams 
17793e3471ccSMark Adams /*@C
1780f1580f4eSBarry Smith   PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is
1781f1580f4eSBarry Smith   called from `PetscFinalize()` automatically.
17823e3471ccSMark Adams 
17833e3471ccSMark Adams   Level: developer
17843e3471ccSMark Adams 
1785*562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscFinalize()`
17863e3471ccSMark Adams @*/
1787d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGFinalizePackage(void)
1788d71ae5a4SJacob Faibussowitsch {
17893e3471ccSMark Adams   PetscFunctionBegin;
17903e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
17919566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&GAMGList));
17923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17933e3471ccSMark Adams }
1794a36cf38bSToby Isaac 
1795a36cf38bSToby Isaac /*@C
1796f1580f4eSBarry Smith   PCGAMGRegister - Register a `PCGAMG` implementation.
1797a36cf38bSToby Isaac 
1798a36cf38bSToby Isaac   Input Parameters:
1799f1580f4eSBarry Smith + type   - string that will be used as the name of the `PCGAMG` type.
1800a36cf38bSToby Isaac - create - function for creating the gamg context.
1801a36cf38bSToby Isaac 
1802f1580f4eSBarry Smith   Level: developer
1803a36cf38bSToby Isaac 
1804*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1805a36cf38bSToby Isaac @*/
1806d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1807d71ae5a4SJacob Faibussowitsch {
1808a36cf38bSToby Isaac   PetscFunctionBegin;
18099566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
18109566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, type, create));
18113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1812a36cf38bSToby Isaac }
18132d776b49SBarry Smith 
18142d776b49SBarry Smith /*@C
18152d776b49SBarry Smith   PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process
18162d776b49SBarry Smith 
18172d776b49SBarry Smith   Input Parameters:
18182d776b49SBarry Smith + pc - the `PCGAMG`
18192d776b49SBarry Smith - A  - the matrix, for any level
18202d776b49SBarry Smith 
18212d776b49SBarry Smith   Output Parameter:
18222d776b49SBarry Smith . G - the graph
18232d776b49SBarry Smith 
18242d776b49SBarry Smith   Level: advanced
18252d776b49SBarry Smith 
1826*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
18272d776b49SBarry Smith @*/
18282d776b49SBarry Smith PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G)
18292d776b49SBarry Smith {
18302d776b49SBarry Smith   PC_MG   *mg      = (PC_MG *)pc->data;
18312d776b49SBarry Smith   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
18322d776b49SBarry Smith 
18332d776b49SBarry Smith   PetscFunctionBegin;
18342d776b49SBarry Smith   PetscCall(pc_gamg->ops->creategraph(pc, A, G));
18353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18362d776b49SBarry Smith }
1837