xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision 8926f930c9b0b9190fd9f60a8cc2d516139e0005)
15b89ad90SMark F. Adams /*
20cd22d39SHong Zhang  GAMG geometric-algebric multigrid PC - Mark Adams 2011
35b89ad90SMark F. Adams  */
4389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h>            /*I "petscpc.h" I*/
518c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/
6f96513f1SMatthew G Knepley 
7c9567895SMark #if defined(PETSC_HAVE_CUDA)
8c9567895SMark   #include <cuda_runtime.h>
9c9567895SMark #endif
10c9567895SMark 
11c9567895SMark #if defined(PETSC_HAVE_HIP)
12c9567895SMark   #include <hip/hip_runtime.h>
13c9567895SMark #endif
14c9567895SMark 
15849bee69SMark Adams PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET];
164555aa8cSStefano Zampini PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
170cbbd2e1SMark F. Adams 
18849bee69SMark Adams // #define GAMG_STAGES
194555aa8cSStefano Zampini #if defined(GAMG_STAGES)
2018c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
21b4fbaa2aSMark F. Adams #endif
22f96513f1SMatthew G Knepley 
230a545947SLisandro Dalcin static PetscFunctionList GAMGList = NULL;
243e3471ccSMark Adams static PetscBool         PCGAMGPackageInitialized;
259d5b6da9SMark F. Adams 
2666976f2fSJacob Faibussowitsch static PetscErrorCode PCReset_GAMG(PC pc)
27d71ae5a4SJacob Faibussowitsch {
28d3d6bff4SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
29d3d6bff4SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
30d3d6bff4SMark F. Adams 
31d3d6bff4SMark F. Adams   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
331c1aac46SBarry Smith   pc_gamg->data_sz = 0;
349566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->orig_data));
355f80ce2aSJacob Faibussowitsch   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) {
3618c3aa7eSMark     mg->min_eigen_DinvA[level] = 0;
3718c3aa7eSMark     mg->max_eigen_DinvA[level] = 0;
3818c3aa7eSMark   }
3918c3aa7eSMark   pc_gamg->emin = 0;
4018c3aa7eSMark   pc_gamg->emax = 0;
41978e3cbaSStefano Zampini   PetscCall(PCReset_MG(pc));
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43a2f3521dSMark F. Adams }
44a2f3521dSMark F. Adams 
455b89ad90SMark F. Adams /*
46c238b0ebSToby Isaac    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
47a147abb0SMark F. Adams      of active processors.
485b89ad90SMark F. Adams 
495b89ad90SMark F. Adams    Input Parameter:
50a2f3521dSMark F. Adams    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
51a2f3521dSMark F. Adams           'pc_gamg->data_sz' are changed via repartitioning/reduction.
529d5b6da9SMark F. Adams    . Amat_fine - matrix on this fine (k) level
53c5bfad50SMark F. Adams    . cr_bs - coarse block size
543530afc2SMark F. Adams    In/Output Parameter:
55a2f3521dSMark F. Adams    . a_P_inout - prolongation operator to the next level (k-->k-1)
56afc97cdcSMark F. Adams    . a_nactive_proc - number of active procs
5711e60469SMark F. Adams    Output Parameter:
583530afc2SMark F. Adams    . a_Amat_crs - coarse matrix that is created (k-1)
595b89ad90SMark F. Adams */
60d71ae5a4SJacob 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)
61d71ae5a4SJacob Faibussowitsch {
629d5b6da9SMark F. Adams   PC_MG      *mg      = (PC_MG *)pc->data;
63486a8d0bSJed Brown   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
64a2f3521dSMark F. Adams   Mat         Cmat, Pold = *a_P_inout;
653b4367a7SBarry Smith   MPI_Comm    comm;
66c5df96a5SBarry Smith   PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc;
673ae0bb68SMark Adams   PetscInt    ncrs_eq, ncrs, f_bs;
685b89ad90SMark F. Adams 
695b89ad90SMark F. Adams   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm));
719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
739566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
74849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
769566063dSJacob Faibussowitsch   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
78849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
79038e3b61SMark F. Adams 
80ce7c7f2fSMark Adams   if (Pcolumnperm) *Pcolumnperm = NULL;
81ce7c7f2fSMark Adams 
823ae0bb68SMark Adams   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
839566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
843ae0bb68SMark Adams   if (pc_gamg->data_cell_rows > 0) {
853ae0bb68SMark Adams     ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows;
8673911c69SBarry Smith   } else {
873ae0bb68SMark Adams     PetscInt bs;
889566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(Cmat, &bs));
893ae0bb68SMark Adams     ncrs = ncrs_eq / bs;
903ae0bb68SMark Adams   }
91c5df96a5SBarry Smith   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
92c9567895SMark   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. */
93c9567895SMark #if defined(PETSC_HAVE_CUDA)
94c9567895SMark     PetscShmComm pshmcomm;
95c9567895SMark     PetscMPIInt  locrank;
96c9567895SMark     MPI_Comm     loccomm;
97c9567895SMark     PetscInt     s_nnodes, r_nnodes, new_new_size;
98c9567895SMark     cudaError_t  cerr;
99c9567895SMark     int          devCount;
1009566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGet(comm, &pshmcomm));
1019566063dSJacob Faibussowitsch     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm));
1029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
103c9567895SMark     s_nnodes = !locrank;
104712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm));
10563a3b9bcSJacob Faibussowitsch     PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes);
106c9567895SMark     devCount = 0;
107c9567895SMark     cerr     = cudaGetDeviceCount(&devCount);
108c9567895SMark     cudaGetLastError();                         /* Reset the last error */
109c9567895SMark     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
110c9567895SMark       new_new_size = r_nnodes * devCount;
111c9567895SMark       new_size     = new_new_size;
11263a3b9bcSJacob 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));
113c9567895SMark     } else {
1149d3446b2SPierre Jolivet       PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.\n", ((PetscObject)pc)->prefix));
115c9567895SMark       goto HEURISTIC;
116c9567895SMark     }
117c9567895SMark #else
118c9567895SMark     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here");
119c9567895SMark #endif
120c9567895SMark   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
1217509f629SMark Adams     if (nactive < pc_gamg->level_reduction_factors[pc_gamg->current_level]) {
1227509f629SMark Adams       new_size = 1;
1237509f629SMark Adams       PetscCall(PetscInfo(pc, "%s: reduction factor too small for %d active processes: reduce to one process\n", ((PetscObject)pc)->prefix, new_size));
1247509f629SMark Adams     } else {
12563a3b9bcSJacob 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]);
126c9567895SMark       new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level];
12763a3b9bcSJacob 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]));
1287509f629SMark Adams     }
129c9567895SMark   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
130c9567895SMark     new_size = 1;
1319566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size));
132c9567895SMark   } else {
133472110cdSMark F. Adams     PetscInt ncrs_eq_glob;
134c9567895SMark #if defined(PETSC_HAVE_CUDA)
135c9567895SMark   HEURISTIC:
136c9567895SMark #endif
1379566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
138a90e85d9SMark Adams     new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
139da81f932SPierre Jolivet     if (!new_size) new_size = 1;                                                       /* not likely, possible? */
140c5df96a5SBarry Smith     else if (new_size >= nactive) new_size = nactive;                                  /* no change, rare */
1419566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size));
142a2f3521dSMark F. Adams   }
1432e3501ffSMark Adams   if (new_size == nactive) {
144ef3f0257SMark Adams     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
145ce7c7f2fSMark Adams     if (new_size < size) {
146ce7c7f2fSMark Adams       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
1479566063dSJacob 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));
148ce7c7f2fSMark Adams       if (pc_gamg->cpu_pin_coarse_grids) {
1499566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE));
1509566063dSJacob Faibussowitsch         PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE));
151ce7c7f2fSMark Adams       }
152ce7c7f2fSMark Adams     }
153ef3f0257SMark Adams     /* we know that the grid structure can be reused in MatPtAP */
1542e3501ffSMark Adams   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
155192c0e8bSMark Adams     PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1;
156885364a3SMark Adams     IS        is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices;
157849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
15871959b99SBarry Smith     nloc_old = ncrs_eq / cr_bs;
15963a3b9bcSJacob 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);
160ce7c7f2fSMark Adams     /* get new_size and rfactor */
161ce7c7f2fSMark Adams     if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
162ce7c7f2fSMark Adams       /* find factor */
163ce7c7f2fSMark Adams       if (new_size == 1) rfactor = size; /* don't modify */
164ce7c7f2fSMark Adams       else {
165ce7c7f2fSMark Adams         PetscReal best_fact = 0.;
166ce7c7f2fSMark Adams         jj                  = -1;
167ce7c7f2fSMark Adams         for (kk = 1; kk <= size; kk++) {
168ce7c7f2fSMark Adams           if (!(size % kk)) { /* a candidate */
169ce7c7f2fSMark Adams             PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size;
170ce7c7f2fSMark Adams             if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */
171ce7c7f2fSMark Adams             if (fact > best_fact) {
1729371c9d4SSatish Balay               best_fact = fact;
1739371c9d4SSatish Balay               jj        = kk;
174ce7c7f2fSMark Adams             }
175ce7c7f2fSMark Adams           }
176ce7c7f2fSMark Adams         }
177ce7c7f2fSMark Adams         if (jj != -1) rfactor = jj;
178ce7c7f2fSMark Adams         else rfactor = 1; /* a prime */
179ce7c7f2fSMark Adams         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
180ce7c7f2fSMark Adams         else expand_factor = rfactor;
181ce7c7f2fSMark Adams       }
182ce7c7f2fSMark Adams       new_size = size / rfactor; /* make new size one that is factor */
1834cdfd227SMark       if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
1844cdfd227SMark         *a_Amat_crs = Cmat;
18563a3b9bcSJacob 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));
186849bee69SMark Adams         PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
1873ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
188ce7c7f2fSMark Adams       }
189ce7c7f2fSMark Adams     }
190a2f3521dSMark F. Adams     /* make 'is_eq_newproc' */
1919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &counts));
192849bee69SMark Adams     if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
1935a9b9e01SMark F. Adams       Mat adj;
194849bee69SMark Adams       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
19563a3b9bcSJacob 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"));
196a2f3521dSMark F. Adams       /* get 'adj' */
197c5bfad50SMark F. Adams       if (cr_bs == 1) {
1989566063dSJacob Faibussowitsch         PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
199806fa848SBarry Smith       } else {
200a2f3521dSMark F. Adams         /* make a scalar matrix to partition (no Stokes here) */
201eb07cef2SMark F. Adams         Mat                tMat;
202a2f3521dSMark F. Adams         PetscInt           Istart_crs, Iend_crs, ncols, jj, Ii;
203b4fbaa2aSMark F. Adams         const PetscScalar *vals;
204b4fbaa2aSMark F. Adams         const PetscInt    *idx;
205a2f3521dSMark F. Adams         PetscInt          *d_nnz, *o_nnz, M, N;
20639d09545SMark Adams         static PetscInt    llev = 0; /* ugly but just used for debugging */
207d9558ea9SBarry Smith         MatType            mtype;
208b4fbaa2aSMark F. Adams 
2099566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz));
2109566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs));
2119566063dSJacob Faibussowitsch         PetscCall(MatGetSize(Cmat, &M, &N));
212c5bfad50SMark F. Adams         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
2139566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL));
214c5bfad50SMark F. Adams           d_nnz[jj] = ncols / cr_bs;
215c5bfad50SMark F. Adams           o_nnz[jj] = ncols / cr_bs;
2169566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL));
2173ae0bb68SMark Adams           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
2183ae0bb68SMark Adams           if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs;
21958471d46SMark F. Adams         }
2206876a03eSMark F. Adams 
2219566063dSJacob Faibussowitsch         PetscCall(MatGetType(Amat_fine, &mtype));
2229566063dSJacob Faibussowitsch         PetscCall(MatCreate(comm, &tMat));
2239566063dSJacob Faibussowitsch         PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE));
2249566063dSJacob Faibussowitsch         PetscCall(MatSetType(tMat, mtype));
2259566063dSJacob Faibussowitsch         PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz));
2269566063dSJacob Faibussowitsch         PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz));
2279566063dSJacob Faibussowitsch         PetscCall(PetscFree2(d_nnz, o_nnz));
228eb07cef2SMark F. Adams 
229a2f3521dSMark F. Adams         for (ii = Istart_crs; ii < Iend_crs; ii++) {
230c5bfad50SMark F. Adams           PetscInt dest_row = ii / cr_bs;
2319566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals));
232eb07cef2SMark F. Adams           for (jj = 0; jj < ncols; jj++) {
233c5bfad50SMark F. Adams             PetscInt    dest_col = idx[jj] / cr_bs;
234eb07cef2SMark F. Adams             PetscScalar v        = 1.0;
2359566063dSJacob Faibussowitsch             PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES));
236eb07cef2SMark F. Adams           }
2379566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals));
238eb07cef2SMark F. Adams         }
2399566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY));
2409566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY));
241eb07cef2SMark F. Adams 
242b4fbaa2aSMark F. Adams         if (llev++ == -1) {
2439371c9d4SSatish Balay           PetscViewer viewer;
2449371c9d4SSatish Balay           char        fname[32];
24563a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev));
2463ba16761SJacob Faibussowitsch           PetscCall(PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer));
2479566063dSJacob Faibussowitsch           PetscCall(MatView(tMat, viewer));
2489566063dSJacob Faibussowitsch           PetscCall(PetscViewerDestroy(&viewer));
249b4fbaa2aSMark F. Adams         }
2509566063dSJacob Faibussowitsch         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
2519566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tMat));
252a2f3521dSMark F. Adams       } /* create 'adj' */
253f150b916SMark F. Adams 
254a2f3521dSMark F. Adams       { /* partition: get newproc_idx */
2555a9b9e01SMark F. Adams         char            prefix[256];
2565a9b9e01SMark F. Adams         const char     *pcpre;
257b4fbaa2aSMark F. Adams         const PetscInt *is_idx;
258b4fbaa2aSMark F. Adams         MatPartitioning mpart;
259a4b7d37bSMark F. Adams         IS              proc_is;
2602f03bc48SMark F. Adams 
2619566063dSJacob Faibussowitsch         PetscCall(MatPartitioningCreate(comm, &mpart));
2629566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
2639566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
2649566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
2659566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
2669566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetFromOptions(mpart));
2679566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, new_size));
2689566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &proc_is));
2699566063dSJacob Faibussowitsch         PetscCall(MatPartitioningDestroy(&mpart));
2705a9b9e01SMark F. Adams 
2715ef31b24SMark F. Adams         /* collect IS info */
2729566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
2739566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(proc_is, &is_idx));
274a2f3521dSMark F. Adams         for (kk = jj = 0; kk < nloc_old; kk++) {
2759371c9d4SSatish Balay           for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ }
2765ef31b24SMark F. Adams         }
2779566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(proc_is, &is_idx));
2789566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&proc_is));
2795ef31b24SMark F. Adams       }
2809566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&adj));
2815a9b9e01SMark F. Adams 
2829566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
2839566063dSJacob Faibussowitsch       PetscCall(PetscFree(newproc_idx));
284849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
28531cb4603SMark Adams     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
286ce7c7f2fSMark Adams       PetscInt targetPE;
28708401ef6SPierre Jolivet       PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen");
28863a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq));
289ce7c7f2fSMark Adams       targetPE = (rank / rfactor) * expand_factor;
2909566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc));
291a2f3521dSMark F. Adams     } /* end simple 'is_eq_newproc' */
292e33ef3b1SMark F. Adams 
29311e60469SMark F. Adams     /*
294a2f3521dSMark F. Adams       Create an index set from the is_eq_newproc index set to indicate the mapping TO
29511e60469SMark F. Adams     */
2969566063dSJacob Faibussowitsch     PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num));
2977700e67bSMark Adams     is_eq_num_prim = is_eq_num;
29811e60469SMark F. Adams     /*
299a2f3521dSMark F. Adams       Determine how many equations/vertices are assigned to each processor
30011e60469SMark F. Adams     */
3019566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(is_eq_newproc, size, counts));
302c5df96a5SBarry Smith     ncrs_eq_new = counts[rank];
3039566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_newproc));
304ce7c7f2fSMark Adams     ncrs_new = ncrs_eq_new / cr_bs;
305a2f3521dSMark F. Adams 
3069566063dSJacob Faibussowitsch     PetscCall(PetscFree(counts));
3076aad120cSJose 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 */
308885364a3SMark Adams     {
309885364a3SMark Adams       Vec             src_crd, dest_crd;
310885364a3SMark 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;
311885364a3SMark Adams       VecScatter      vecscat;
312885364a3SMark Adams       PetscScalar    *array;
313885364a3SMark Adams       IS              isscat;
314a2f3521dSMark F. Adams       /* move data (for primal equations only) */
31522063be5SMark F. Adams       /* Create a vector to contain the newly ordered element information */
3169566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &dest_crd));
3179566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE));
3189566063dSJacob Faibussowitsch       PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */
31911e60469SMark F. Adams       /*
3209d5b6da9SMark F. Adams         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
321c5bfad50SMark F. Adams         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
32211e60469SMark F. Adams       */
3239566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx));
3249566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(is_eq_num_prim, &idx));
3253ae0bb68SMark Adams       for (ii = 0, jj = 0; ii < ncrs; ii++) {
326c5bfad50SMark F. Adams         PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */
327a2f3521dSMark F. Adams         for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk;
32811e60469SMark F. Adams       }
3299566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(is_eq_num_prim, &idx));
3309566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat));
3319566063dSJacob Faibussowitsch       PetscCall(PetscFree(tidx));
33211e60469SMark F. Adams       /*
33311e60469SMark F. Adams         Create a vector to contain the original vertex information for each element
33411e60469SMark F. Adams       */
3359566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd));
3369d5b6da9SMark F. Adams       for (jj = 0; jj < ndata_cols; jj++) {
3373ae0bb68SMark Adams         const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows;
3383ae0bb68SMark Adams         for (ii = 0; ii < ncrs; ii++) {
3399d5b6da9SMark F. Adams           for (kk = 0; kk < ndata_rows; kk++) {
340a2f3521dSMark F. Adams             PetscInt    ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj;
341c8b0795cSMark F. Adams             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
3429566063dSJacob Faibussowitsch             PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES));
343d3d6bff4SMark F. Adams           }
344038e3b61SMark F. Adams         }
345eb07cef2SMark F. Adams       }
3469566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(src_crd));
3479566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(src_crd));
34811e60469SMark F. Adams       /*
34911e60469SMark F. Adams         Scatter the element vertex information (still in the original vertex ordering)
35011e60469SMark F. Adams         to the correct processor
35111e60469SMark F. Adams       */
3529566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat));
3539566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isscat));
3549566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD));
3559566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD));
3569566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&vecscat));
3579566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&src_crd));
35811e60469SMark F. Adams       /*
35911e60469SMark F. Adams         Put the element vertex data into a new allocation of the gdata->ele
36011e60469SMark F. Adams       */
3619566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
3629566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data));
3632fa5cd67SKarl Rupp 
3643ae0bb68SMark Adams       pc_gamg->data_sz = node_data_sz * ncrs_new;
3653ae0bb68SMark Adams       strideNew        = ncrs_new * ndata_rows;
3662fa5cd67SKarl Rupp 
3679566063dSJacob Faibussowitsch       PetscCall(VecGetArray(dest_crd, &array));
3689d5b6da9SMark F. Adams       for (jj = 0; jj < ndata_cols; jj++) {
3693ae0bb68SMark Adams         for (ii = 0; ii < ncrs_new; ii++) {
3709d5b6da9SMark F. Adams           for (kk = 0; kk < ndata_rows; kk++) {
371a2f3521dSMark F. Adams             PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj;
372c8b0795cSMark F. Adams             pc_gamg->data[ix] = PetscRealPart(array[jx]);
373d3d6bff4SMark F. Adams           }
374038e3b61SMark F. Adams         }
375038e3b61SMark F. Adams       }
3769566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(dest_crd, &array));
3779566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&dest_crd));
378885364a3SMark Adams     }
379a2f3521dSMark F. Adams     /* move A and P (columns) with new layout */
38011e60469SMark F. Adams     /*
3817dae84e0SHong Zhang       Invert for MatCreateSubMatrix
38211e60469SMark F. Adams     */
3839566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
3849566063dSJacob Faibussowitsch     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
3859566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
3869371c9d4SSatish Balay     if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ }
3873cb8563fSToby Isaac     if (Pcolumnperm) {
3889566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
3893cb8563fSToby Isaac       *Pcolumnperm = new_eq_indices;
3903cb8563fSToby Isaac     }
3919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_eq_num));
392849bee69SMark Adams 
393a2f3521dSMark F. Adams     /* 'a_Amat_crs' output */
394a2f3521dSMark F. Adams     {
395a2f3521dSMark F. Adams       Mat       mat;
396b94d7dedSBarry Smith       PetscBool isset, isspd, isher;
39790db8557SMark Adams #if !defined(PETSC_USE_COMPLEX)
398b94d7dedSBarry Smith       PetscBool issym;
399b94d7dedSBarry Smith #endif
400b94d7dedSBarry Smith 
401b94d7dedSBarry Smith       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
402b94d7dedSBarry Smith       PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ?
403b94d7dedSBarry Smith       if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd));
404b94d7dedSBarry Smith       else {
405b94d7dedSBarry Smith         PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher));
406b94d7dedSBarry Smith         if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher));
407b94d7dedSBarry Smith         else {
408b94d7dedSBarry Smith #if !defined(PETSC_USE_COMPLEX)
409b94d7dedSBarry Smith           PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym));
410b94d7dedSBarry Smith           if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym));
41190db8557SMark Adams #endif
41290db8557SMark Adams         }
41390db8557SMark Adams       }
414a2f3521dSMark F. Adams       *a_Amat_crs = mat;
415a2f3521dSMark F. Adams     }
4169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cmat));
417a2f3521dSMark F. Adams 
41811e60469SMark F. Adams     /* prolongator */
41911e60469SMark F. Adams     {
42011e60469SMark F. Adams       IS       findices;
421a2f3521dSMark F. Adams       PetscInt Istart, Iend;
422a2f3521dSMark F. Adams       Mat      Pnew;
42362294041SBarry Smith 
4249566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
4259566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices));
4269566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(findices, f_bs));
4279566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
4289566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&findices));
4299566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
430c5bfad50SMark F. Adams 
4319566063dSJacob Faibussowitsch       PetscCall(MatDestroy(a_P_inout));
432a2f3521dSMark F. Adams 
433a2f3521dSMark F. Adams       /* output - repartitioned */
434a2f3521dSMark F. Adams       *a_P_inout = Pnew;
435e33ef3b1SMark F. Adams     }
4369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&new_eq_indices));
4375b89ad90SMark F. Adams 
438c5df96a5SBarry Smith     *a_nactive_proc = new_size; /* output */
439ce7c7f2fSMark Adams 
440ce7c7f2fSMark Adams     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
441ce7c7f2fSMark Adams     if (pc_gamg->cpu_pin_coarse_grids) {
442ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
4438bca76a6SMark Adams       static PetscInt llev = 2;
44463a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Pinning level %" PetscInt_FMT " to the CPU\n", ((PetscObject)pc)->prefix, llev++));
445ce7c7f2fSMark Adams #endif
4469566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE));
4479566063dSJacob Faibussowitsch       PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE));
448adf5291fSStefano Zampini       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
449ce7c7f2fSMark Adams         Mat         A = *a_Amat_crs, P = *a_P_inout;
450ce7c7f2fSMark Adams         PetscMPIInt size;
4519566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
452ce7c7f2fSMark Adams         if (size > 1) {
453ce7c7f2fSMark Adams           Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *p = (Mat_MPIAIJ *)P->data;
4549566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(a->lvec, PETSC_TRUE));
4559566063dSJacob Faibussowitsch           PetscCall(VecBindToCPU(p->lvec, PETSC_TRUE));
456ce7c7f2fSMark Adams         }
457ce7c7f2fSMark Adams       }
458ce7c7f2fSMark Adams     }
459849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
460849bee69SMark Adams   } // processor reduce
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4625b89ad90SMark F. Adams }
4635b89ad90SMark F. Adams 
464bae903cbSmarkadams4 // used in GEO
465d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2)
466d71ae5a4SJacob Faibussowitsch {
4674b1575e2SStefano Zampini   const char *prefix;
4684b1575e2SStefano Zampini   char        addp[32];
4694b1575e2SStefano Zampini   PC_MG      *mg      = (PC_MG *)a_pc->data;
4704b1575e2SStefano Zampini   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
4714b1575e2SStefano Zampini 
4724b1575e2SStefano Zampini   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(a_pc, &prefix));
47463a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(a_pc, "%s: Square Graph on level %" PetscInt_FMT "\n", ((PetscObject)a_pc)->prefix, pc_gamg->current_level + 1));
4759566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(Gmat1, Gmat1, NULL, Gmat2));
4769566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(*Gmat2, prefix));
47763a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_square_%" PetscInt_FMT "_", pc_gamg->current_level));
4789566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(*Gmat2, addp));
479b94d7dedSBarry Smith   if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) {
4809566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AB));
481b4da3a1bSStefano Zampini   } else {
4829566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Gmat1, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
4839566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(*Gmat2, MATPRODUCT_AtB));
484b4da3a1bSStefano Zampini   }
4859566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(*Gmat2));
4869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0));
4879566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(*Gmat2));
4889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0], 0, 0, 0, 0));
4899566063dSJacob Faibussowitsch   PetscCall(MatProductClear(*Gmat2));
4904b1575e2SStefano Zampini   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
4914b1575e2SStefano Zampini   (*Gmat2)->assembled = PETSC_TRUE;
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4934b1575e2SStefano Zampini }
4944b1575e2SStefano Zampini 
4955b89ad90SMark F. Adams /*
4965b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
4975b89ad90SMark F. Adams                     by setting data structures and options.
4985b89ad90SMark F. Adams 
4995b89ad90SMark F. Adams    Input Parameter:
5005b89ad90SMark F. Adams .  pc - the preconditioner context
5015b89ad90SMark F. Adams 
5025b89ad90SMark F. Adams */
50366976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_GAMG(PC pc)
504d71ae5a4SJacob Faibussowitsch {
5059d5b6da9SMark F. Adams   PC_MG      *mg      = (PC_MG *)pc->data;
5065b89ad90SMark F. Adams   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
5072adcac29SMark F. Adams   Mat         Pmat    = pc->pmat;
50818c3aa7eSMark   PetscInt    fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS];
5093b4367a7SBarry Smith   MPI_Comm    comm;
510c5df96a5SBarry Smith   PetscMPIInt rank, size, nactivepe;
51118c3aa7eSMark   Mat         Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS];
51218c3aa7eSMark   IS         *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
5134279555eSSatish Balay   PetscBool   is_last = PETSC_FALSE;
5144279555eSSatish Balay #if defined(PETSC_USE_INFO)
515a2f3521dSMark F. Adams   PetscLogDouble nnz0 = 0., nnztot = 0.;
516569f4572SMark Adams   MatInfo        info;
5174279555eSSatish Balay #endif
5185ef31b24SMark F. Adams 
5195b89ad90SMark F. Adams   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
5219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
523849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
5248abdc6daSStefano Zampini   if (pc->setupcalled) {
5258abdc6daSStefano Zampini     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
526878e152fSMark F. Adams       /* reset everything */
5279566063dSJacob Faibussowitsch       PetscCall(PCReset_MG(pc));
528878e152fSMark F. Adams       pc->setupcalled = 0;
529806fa848SBarry Smith     } else {
53084d3f75bSMark F. Adams       PC_MG_Levels **mglevels = mg->levels;
53103a628feSMark F. Adams       /* just do Galerkin grids */
53258471d46SMark F. Adams       Mat B, dA, dB;
5339d5b6da9SMark F. Adams       if (pc_gamg->Nlevels > 1) {
5344555aa8cSStefano Zampini         PetscInt gl;
53558471d46SMark F. Adams         /* currently only handle case where mat and pmat are the same on coarser levels */
5369566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB));
53758471d46SMark F. Adams         /* (re)set to get dirty flag */
5389566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB));
53958471d46SMark F. Adams 
5404555aa8cSStefano Zampini         for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) {
5418abdc6daSStefano Zampini           MatReuse reuse = MAT_INITIAL_MATRIX;
542849bee69SMark Adams #if defined(GAMG_STAGES)
543849bee69SMark Adams           PetscCall(PetscLogStagePush(gamg_stages[gl]));
544849bee69SMark Adams #endif
5458abdc6daSStefano Zampini           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
5469566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(mglevels[level]->smoothd, NULL, &B));
5478abdc6daSStefano Zampini           if (B->product) {
548ad540459SPierre Jolivet             if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX;
5498abdc6daSStefano Zampini           }
5509566063dSJacob Faibussowitsch           if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A));
5518abdc6daSStefano Zampini           if (reuse == MAT_REUSE_MATRIX) {
552e02fb3cdSMark Adams             PetscCall(PetscInfo(pc, "%s: RAP after initial setup, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
5538abdc6daSStefano Zampini           } else {
554e02fb3cdSMark Adams             PetscCall(PetscInfo(pc, "%s: RAP after initial setup, with repartitioning (new matrix) level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
5558abdc6daSStefano Zampini           }
5569566063dSJacob Faibussowitsch           PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0));
5579566063dSJacob Faibussowitsch           PetscCall(MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B));
5589566063dSJacob Faibussowitsch           PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0));
55963b77682SMark Adams           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
5609566063dSJacob Faibussowitsch           PetscCall(KSPSetOperators(mglevels[level]->smoothd, B, B));
561e1cf1444SMark Adams           // check for redoing eigen estimates
562e1cf1444SMark Adams           if (pc_gamg->recompute_esteig) {
563e1cf1444SMark Adams             PetscBool ischeb;
564e1cf1444SMark Adams             KSP       smoother;
565e1cf1444SMark Adams             PetscCall(PCMGGetSmoother(pc, level + 1, &smoother));
566e1cf1444SMark Adams             PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
567e1cf1444SMark Adams             if (ischeb) {
568e1cf1444SMark Adams               KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
569e1cf1444SMark Adams               cheb->emin_provided = 0;
570e1cf1444SMark Adams               cheb->emax_provided = 0;
571e1cf1444SMark Adams             }
572e1cf1444SMark Adams             /* we could call PetscCall(KSPChebyshevSetEigenvalues(smoother, 0, 0)); but the logic does not work properly */
573e1cf1444SMark Adams           }
574e1cf1444SMark Adams           // inc
57558471d46SMark F. Adams           dB = B;
576849bee69SMark Adams #if defined(GAMG_STAGES)
577849bee69SMark Adams           PetscCall(PetscLogStagePop());
578849bee69SMark Adams #endif
57958471d46SMark F. Adams         }
5805f8cf99dSMark F. Adams       }
581d5280255SMark F. Adams 
5829566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
583849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
5843ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
585eb07cef2SMark F. Adams     }
586878e152fSMark F. Adams   }
587f6536408SMark F. Adams 
588878e152fSMark F. Adams   if (!pc_gamg->data) {
589878e152fSMark F. Adams     if (pc_gamg->orig_data) {
5909566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(Pmat, &bs));
5919566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
5922fa5cd67SKarl Rupp 
593878e152fSMark F. Adams       pc_gamg->data_sz        = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols;
594878e152fSMark F. Adams       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
595878e152fSMark F. Adams       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
5962fa5cd67SKarl Rupp 
5979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
598878e152fSMark F. Adams       for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
599806fa848SBarry Smith     } else {
6005f80ce2aSJacob Faibussowitsch       PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data");
6019566063dSJacob Faibussowitsch       PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat));
6029d5b6da9SMark F. Adams     }
603878e152fSMark F. Adams   }
604878e152fSMark F. Adams 
605878e152fSMark F. Adams   /* cache original data for reuse */
6061c1aac46SBarry Smith   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
6079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
608878e152fSMark F. Adams     for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
609878e152fSMark F. Adams     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
610878e152fSMark F. Adams     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
611878e152fSMark F. Adams   }
612038e3b61SMark F. Adams 
613302f38e8SMark F. Adams   /* get basic dims */
6149566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Pmat, &bs));
6159566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Pmat, &M, &N));
61684d3f75bSMark F. Adams 
6174279555eSSatish Balay #if defined(PETSC_USE_INFO)
6189566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */
619569f4572SMark Adams   nnz0   = info.nz_used;
620569f4572SMark Adams   nnztot = info.nz_used;
6214279555eSSatish Balay #endif
622bae903cbSmarkadams4   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));
623569f4572SMark Adams 
624a2f3521dSMark F. Adams   /* Get A_i and R_i */
62562294041SBarry Smith   for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) {
6269ab59c8bSMark Adams     pc_gamg->current_level = level;
6279c15c1aeSMark Adams     PetscCheck(level < PETSC_MG_MAXLEVELS - 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level + 1);
6285b89ad90SMark F. Adams     level1 = level + 1;
6294555aa8cSStefano Zampini #if defined(GAMG_STAGES)
630849bee69SMark Adams     if (!gamg_stages[level]) {
631849bee69SMark Adams       char str[32];
632a364092eSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level));
633849bee69SMark Adams       PetscCall(PetscLogStageRegister(str, &gamg_stages[level]));
634849bee69SMark Adams     }
6359566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(gamg_stages[level]));
636b4fbaa2aSMark F. Adams #endif
637c8b0795cSMark F. Adams     { /* construct prolongator */
638*8926f930SMark Adams       Mat               Gmat, mat;
6390cbbd2e1SMark F. Adams       PetscCoarsenData *agg_lists;
6407700e67bSMark Adams       Mat               Prol11;
641c8b0795cSMark F. Adams 
6422d776b49SBarry Smith       PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat));
6439c15c1aeSMark Adams       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); // Gmat may have ghosts for QR aggregates not in matrix
644*8926f930SMark Adams       PetscCall(PetscCDGetMat(agg_lists, &mat));
645*8926f930SMark Adams       if (!mat) PetscCall(PetscCDSetMat(agg_lists, Gmat));
646*8926f930SMark Adams       PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], agg_lists, &Prol11));
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 */
653*8926f930SMark Adams         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); // column size
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) {
661*8926f930SMark Adams           PetscInt bs;
662*8926f930SMark Adams           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // row block size
6639c15c1aeSMark Adams           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
6649c15c1aeSMark Adams           PetscCall(PetscInfo(pc, "%d: %" PetscInt_FMT " ASM local domains,  bs = %d\n", (int)level, nASMBlocksArr[level], (int)bs));
665*8926f930SMark Adams         } else if (pc_gamg->asm_hem_aggs) {
666*8926f930SMark Adams           MatCoarsen  crs;
667*8926f930SMark Adams           const char *prefix;
668*8926f930SMark Adams           PetscInt    bs;
669*8926f930SMark Adams           PetscCall(PetscCDGetMat(agg_lists, &mat));
670*8926f930SMark Adams           if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck)
671*8926f930SMark Adams           PetscCall(PetscCDDestroy(agg_lists));
672*8926f930SMark Adams           PetscCall(PetscInfo(pc, "HEM ASM passes = %d\n", (int)pc_gamg->asm_hem_aggs));
673*8926f930SMark Adams           PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &crs));
674*8926f930SMark Adams           PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
675*8926f930SMark Adams           PetscCall(PetscObjectSetOptionsPrefix((PetscObject)crs, prefix));
676*8926f930SMark Adams           PetscCall(MatCoarsenSetFromOptions(crs)); // get strength args
677*8926f930SMark Adams           PetscCall(MatCoarsenSetType(crs, MATCOARSENHEM));
678*8926f930SMark Adams           PetscCall(MatCoarsenSetMaximumIterations(crs, pc_gamg->asm_hem_aggs));
679*8926f930SMark Adams           PetscCall(MatCoarsenSetAdjacency(crs, Gmat));
680*8926f930SMark Adams           PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
681*8926f930SMark Adams           PetscCall(MatCoarsenApply(crs));
682*8926f930SMark Adams           PetscCall(MatCoarsenViewFromOptions(crs, NULL, "-agg_hem_mat_coarsen_view"));
683*8926f930SMark Adams           PetscCall(MatCoarsenGetData(crs, &agg_lists)); /* output */
684*8926f930SMark Adams           PetscCall(MatCoarsenDestroy(&crs));
685*8926f930SMark Adams           // create aggregates
686*8926f930SMark Adams           PetscCall(MatGetBlockSizes(Aarr[level], &bs, NULL)); // row block size
687*8926f930SMark Adams           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
688ffc955d6SMark F. Adams         }
6899566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
6909566063dSJacob Faibussowitsch         PetscCall(MatSetOptionsPrefix(Prol11, prefix));
6919566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level));
6929566063dSJacob Faibussowitsch         PetscCall(MatAppendOptionsPrefix(Prol11, addp));
69391f31d3dSStefano Zampini         /* Always generate the transpose with CUDA
694f7df55f0SStefano Zampini            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
6959566063dSJacob Faibussowitsch         PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
6969566063dSJacob Faibussowitsch         PetscCall(MatSetFromOptions(Prol11));
6974bde40a0SMark Adams         Parr[level1] = Prol11;
6984bde40a0SMark Adams       } else Parr[level1] = NULL; /* failed to coarsen */
699*8926f930SMark Adams       PetscCall(PetscCDGetMat(agg_lists, &mat));
700*8926f930SMark Adams       if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck)
7019566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat));
7029566063dSJacob Faibussowitsch       PetscCall(PetscCDDestroy(agg_lists));
703a2f3521dSMark F. Adams     }                           /* construct prolongator scope */
7047f66b68fSMark Adams     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
705171cca9aSMark Adams     if (!Parr[level1]) {        /* failed to coarsen */
70663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
7074555aa8cSStefano Zampini #if defined(GAMG_STAGES)
7089566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePop());
709a90e85d9SMark Adams #endif
710c8b0795cSMark F. Adams       break;
711c8b0795cSMark F. Adams     }
7129566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
7132472a847SBarry Smith     PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?");
714171cca9aSMark Adams     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
7150e2909e1SMark Adams     if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE;
716849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
7179566063dSJacob Faibussowitsch     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
718849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
719a2f3521dSMark F. Adams 
7209566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
7214279555eSSatish Balay #if defined(PETSC_USE_INFO)
7229566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
723569f4572SMark Adams     nnztot += info.nz_used;
7244279555eSSatish Balay #endif
725bae903cbSmarkadams4     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));
726569f4572SMark Adams 
7274555aa8cSStefano Zampini #if defined(GAMG_STAGES)
7289566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
729b4fbaa2aSMark F. Adams #endif
730a90e85d9SMark Adams     /* stop if one node or one proc -- could pull back for singular problems */
7319ab59c8bSMark Adams     if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) {
73263a3b9bcSJacob 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));
733a90e85d9SMark Adams       level++;
734a90e85d9SMark Adams       break;
735a90e85d9SMark Adams     }
736c8b0795cSMark F. Adams   } /* levels */
7379566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->data));
738c8b0795cSMark F. Adams 
739ba4ceb06Smarkadams4   PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0));
7409d5b6da9SMark F. Adams   pc_gamg->Nlevels = level + 1;
7415b89ad90SMark F. Adams   fine_level       = level;
7429566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL));
7435b89ad90SMark F. Adams 
74462294041SBarry Smith   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
7450ed2132dSStefano Zampini 
746d5280255SMark F. Adams     /* set default smoothers & set operators */
74762294041SBarry Smith     for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) {
748ffc955d6SMark F. Adams       KSP smoother;
749ffc955d6SMark F. Adams       PC  subpc;
750a2f3521dSMark F. Adams 
7519566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7529566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(smoother, &subpc));
753ffc955d6SMark F. Adams 
7549566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
755a2f3521dSMark F. Adams       /* set ops */
7569566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
7579566063dSJacob Faibussowitsch       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1]));
758a2f3521dSMark F. Adams 
759a2f3521dSMark F. Adams       /* set defaults */
7609566063dSJacob Faibussowitsch       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
761a2f3521dSMark F. Adams 
7620c3bc534SBarry Smith       /* set blocks for ASM smoother that uses the 'aggregates' */
763*8926f930SMark Adams       if (pc_gamg->use_aggs_in_asm || pc_gamg->asm_hem_aggs) {
7642d3561bbSSatish Balay         PetscInt sz;
7657a28f3e5SMark Adams         IS      *iss;
766a2f3521dSMark F. Adams 
7672d3561bbSSatish Balay         sz  = nASMBlocksArr[level];
7687a28f3e5SMark Adams         iss = ASMLocalIDsArr[level];
7699566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCASM));
7709566063dSJacob Faibussowitsch         PetscCall(PCASMSetOverlap(subpc, 0));
7719566063dSJacob Faibussowitsch         PetscCall(PCASMSetType(subpc, PC_ASM_BASIC));
7727f66b68fSMark Adams         if (!sz) {
773ffc955d6SMark F. Adams           IS is;
7749566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
7759566063dSJacob Faibussowitsch           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
7769566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
777806fa848SBarry Smith         } else {
778a94c3b12SMark F. Adams           PetscInt kk;
7792eab5cd7Smarkadams4           PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL));
78048a46eb9SPierre Jolivet           for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk]));
7819566063dSJacob Faibussowitsch           PetscCall(PetscFree(iss));
782ffc955d6SMark F. Adams         }
7830298fd71SBarry Smith         ASMLocalIDsArr[level] = NULL;
784ffc955d6SMark F. Adams         nASMBlocksArr[level]  = 0;
785806fa848SBarry Smith       } else {
7869566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCJACOBI));
787ffc955d6SMark F. Adams       }
788d5280255SMark F. Adams     }
789d5280255SMark F. Adams     {
790d5280255SMark F. Adams       /* coarse grid */
7919371c9d4SSatish Balay       KSP      smoother, *k2;
7929371c9d4SSatish Balay       PC       subpc, pc2;
7939371c9d4SSatish Balay       PetscInt ii, first;
7949371c9d4SSatish Balay       Mat      Lmat = Aarr[(level = pc_gamg->Nlevels - 1)];
7959371c9d4SSatish Balay       lidx          = 0;
7960ed2132dSStefano Zampini 
7979566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
7989566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
799cf8ae1d3SMark Adams       if (!pc_gamg->use_parallel_coarse_grid_solver) {
8009566063dSJacob Faibussowitsch         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
8019566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(smoother, &subpc));
8029566063dSJacob Faibussowitsch         PetscCall(PCSetType(subpc, PCBJACOBI));
8039566063dSJacob Faibussowitsch         PetscCall(PCSetUp(subpc));
8049566063dSJacob Faibussowitsch         PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2));
80563a3b9bcSJacob Faibussowitsch         PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii);
8069566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(k2[0], &pc2));
8079566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc2, PCLU));
8089566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS));
8099566063dSJacob Faibussowitsch         PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1));
8109566063dSJacob Faibussowitsch         PetscCall(KSPSetType(k2[0], KSPPREONLY));
811d5280255SMark F. Adams       }
812cf8ae1d3SMark Adams     }
813d5280255SMark F. Adams 
814d5280255SMark F. Adams     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
815d0609cedSBarry Smith     PetscObjectOptionsBegin((PetscObject)pc);
816dbbe0bcdSBarry Smith     PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject));
817d0609cedSBarry Smith     PetscOptionsEnd();
8189566063dSJacob Faibussowitsch     PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
819d5280255SMark F. Adams 
820f1580f4eSBarry Smith     /* set cheby eigen estimates from SA to use in the solver */
8217e6512fdSJed Brown     if (pc_gamg->use_sa_esteig) {
82218c3aa7eSMark       for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) {
82318c3aa7eSMark         KSP       smoother;
82418c3aa7eSMark         PetscBool ischeb;
8250ed2132dSStefano Zampini 
8269566063dSJacob Faibussowitsch         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
8279566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
82818c3aa7eSMark         if (ischeb) {
82918c3aa7eSMark           KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
8300ed2132dSStefano Zampini 
8312de708cbSJed Brown           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
832efe053fcSJed Brown           if (mg->max_eigen_DinvA[level] > 0) {
8337e6512fdSJed Brown             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
8347e6512fdSJed Brown             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
83518c3aa7eSMark             PetscReal emax, emin;
8360ed2132dSStefano Zampini 
83718c3aa7eSMark             emin = mg->min_eigen_DinvA[level];
83818c3aa7eSMark             emax = mg->max_eigen_DinvA[level];
83963a3b9bcSJacob 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));
8400a94a983SJed Brown             cheb->emin_provided = emin;
8410a94a983SJed Brown             cheb->emax_provided = emax;
84218c3aa7eSMark           }
84318c3aa7eSMark         }
84418c3aa7eSMark       }
84518c3aa7eSMark     }
8460ed2132dSStefano Zampini 
8479566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8480ed2132dSStefano Zampini 
849d5280255SMark F. Adams     /* clean up */
850d5280255SMark F. Adams     for (level = 1; level < pc_gamg->Nlevels; level++) {
8519566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Parr[level]));
8529566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Aarr[level]));
8535b89ad90SMark F. Adams     }
854806fa848SBarry Smith   } else {
8555f8cf99dSMark F. Adams     KSP smoother;
8560ed2132dSStefano Zampini 
85763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix));
8589566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
8599566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
8609566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPPREONLY));
8619566063dSJacob Faibussowitsch     PetscCall(PCSetUp_MG(pc));
8625f8cf99dSMark F. Adams   }
863849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
8643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8655b89ad90SMark F. Adams }
8665b89ad90SMark F. Adams 
8675b89ad90SMark F. Adams /*
8685b89ad90SMark F. Adams  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
8695b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
8705b89ad90SMark F. Adams 
8715b89ad90SMark F. Adams    Input Parameter:
8725b89ad90SMark F. Adams .  pc - the preconditioner context
8735b89ad90SMark F. Adams 
8745b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
8755b89ad90SMark F. Adams */
876d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_GAMG(PC pc)
877d71ae5a4SJacob Faibussowitsch {
8785b89ad90SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
8795b89ad90SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
8805b89ad90SMark F. Adams 
8815b89ad90SMark F. Adams   PetscFunctionBegin;
8829566063dSJacob Faibussowitsch   PetscCall(PCReset_GAMG(pc));
8831baa6e33SBarry Smith   if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc));
8849566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->ops));
8859566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
8869566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg));
8872e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL));
8882e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL));
8892e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL));
8902e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL));
8912e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL));
892e1cf1444SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", NULL));
8932e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL));
8942e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL));
895d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", NULL));
8962e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL));
8972e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL));
8982e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL));
8992e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL));
9002e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL));
9012e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL));
9022e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL));
9032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL));
904*8926f930SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", NULL));
9059566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
9063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9075b89ad90SMark F. Adams }
9085b89ad90SMark F. Adams 
909676e1743SMark F. Adams /*@
910f1580f4eSBarry Smith   PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG`
911676e1743SMark F. Adams 
912c3339decSBarry Smith   Logically Collective
913676e1743SMark F. Adams 
914676e1743SMark F. Adams   Input Parameters:
9151cc46a46SBarry Smith + pc - the preconditioner context
9161cc46a46SBarry Smith - n  - the number of equations
917676e1743SMark F. Adams 
918676e1743SMark F. Adams   Options Database Key:
919147403d9SBarry Smith . -pc_gamg_process_eq_limit <limit> - set the limit
920676e1743SMark F. Adams 
92120f4b53cSBarry Smith   Level: intermediate
92220f4b53cSBarry Smith 
923f1580f4eSBarry Smith   Note:
924f1580f4eSBarry 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
925cab9ed1eSBarry Smith   that has degrees of freedom
926cab9ed1eSBarry Smith 
927562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
928676e1743SMark F. Adams @*/
929d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
930d71ae5a4SJacob Faibussowitsch {
931676e1743SMark F. Adams   PetscFunctionBegin;
932676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
933cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n));
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935676e1743SMark F. Adams }
936676e1743SMark F. Adams 
937d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
938d71ae5a4SJacob Faibussowitsch {
939c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
940c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
941676e1743SMark F. Adams 
942676e1743SMark F. Adams   PetscFunctionBegin;
9439d5b6da9SMark F. Adams   if (n > 0) pc_gamg->min_eq_proc = n;
9443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
945676e1743SMark F. Adams }
946676e1743SMark F. Adams 
947389730f3SMark F. Adams /*@
948f1580f4eSBarry Smith   PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG`
949389730f3SMark F. Adams 
950c3339decSBarry Smith   Collective
951389730f3SMark F. Adams 
952389730f3SMark F. Adams   Input Parameters:
9531cc46a46SBarry Smith + pc - the preconditioner context
9541cc46a46SBarry Smith - n  - maximum number of equations to aim for
955389730f3SMark F. Adams 
956389730f3SMark F. Adams   Options Database Key:
957147403d9SBarry Smith . -pc_gamg_coarse_eq_limit <limit> - set the limit
958389730f3SMark F. Adams 
95920f4b53cSBarry Smith   Level: intermediate
96020f4b53cSBarry Smith 
961f1580f4eSBarry Smith   Note:
962f1580f4eSBarry Smith   For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
96374329af1SBarry Smith   has less than 1000 unknowns.
96474329af1SBarry Smith 
965562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
966389730f3SMark F. Adams @*/
967d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
968d71ae5a4SJacob Faibussowitsch {
969389730f3SMark F. Adams   PetscFunctionBegin;
970389730f3SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
971cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n));
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
973389730f3SMark F. Adams }
974389730f3SMark F. Adams 
975d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
976d71ae5a4SJacob Faibussowitsch {
977389730f3SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
978389730f3SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
979389730f3SMark F. Adams 
980389730f3SMark F. Adams   PetscFunctionBegin;
9819d5b6da9SMark F. Adams   if (n > 0) pc_gamg->coarse_eq_limit = n;
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
983389730f3SMark F. Adams }
984389730f3SMark F. Adams 
985676e1743SMark F. Adams /*@
986f1580f4eSBarry Smith   PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use
987676e1743SMark F. Adams 
988c3339decSBarry Smith   Collective
989676e1743SMark F. Adams 
990676e1743SMark F. Adams   Input Parameters:
9911cc46a46SBarry Smith + pc - the preconditioner context
992f1580f4eSBarry Smith - n  - `PETSC_TRUE` or `PETSC_FALSE`
993676e1743SMark F. Adams 
994676e1743SMark F. Adams   Options Database Key:
995147403d9SBarry Smith . -pc_gamg_repartition <true,false> - turn on the repartitioning
996676e1743SMark F. Adams 
99720f4b53cSBarry Smith   Level: intermediate
99820f4b53cSBarry Smith 
999f1580f4eSBarry Smith   Note:
1000f1580f4eSBarry Smith   This will generally improve the loading balancing of the work on each level
1001cab9ed1eSBarry Smith 
1002562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
1003676e1743SMark F. Adams @*/
1004d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
1005d71ae5a4SJacob Faibussowitsch {
1006676e1743SMark F. Adams   PetscFunctionBegin;
1007676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1008cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n));
10093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1010676e1743SMark F. Adams }
1011676e1743SMark F. Adams 
1012d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
1013d71ae5a4SJacob Faibussowitsch {
1014c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1015c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1016676e1743SMark F. Adams 
1017676e1743SMark F. Adams   PetscFunctionBegin;
10189d5b6da9SMark F. Adams   pc_gamg->repart = n;
10193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1020676e1743SMark F. Adams }
1021676e1743SMark F. Adams 
1022dfd5c07aSMark F. Adams /*@
1023f1580f4eSBarry Smith   PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process
102418c3aa7eSMark 
1025c3339decSBarry Smith   Collective
102618c3aa7eSMark 
102718c3aa7eSMark   Input Parameters:
102818c3aa7eSMark + pc - the preconditioner context
1029e1cf1444SMark Adams - b  - flag
103018c3aa7eSMark 
103118c3aa7eSMark   Options Database Key:
1032147403d9SBarry Smith . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
103318c3aa7eSMark 
103420f4b53cSBarry Smith   Level: advanced
103520f4b53cSBarry Smith 
103618c3aa7eSMark   Notes:
10377e6512fdSJed 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$.
1038f1580f4eSBarry Smith   Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
1039f1580f4eSBarry Smith   If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process
1040f1580f4eSBarry Smith   This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used.
1041efe053fcSJed Brown   It became default in PETSc 3.17.
104218c3aa7eSMark 
1043562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `PCGAMGSetRecomputeEstEig()`
104418c3aa7eSMark @*/
1045e1cf1444SMark Adams PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool b)
1046d71ae5a4SJacob Faibussowitsch {
104718c3aa7eSMark   PetscFunctionBegin;
104818c3aa7eSMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1049e1cf1444SMark Adams   PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, b));
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105118c3aa7eSMark }
105218c3aa7eSMark 
1053e1cf1444SMark Adams static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool b)
1054d71ae5a4SJacob Faibussowitsch {
105518c3aa7eSMark   PC_MG   *mg      = (PC_MG *)pc->data;
105618c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
105718c3aa7eSMark 
105818c3aa7eSMark   PetscFunctionBegin;
1059e1cf1444SMark Adams   pc_gamg->use_sa_esteig = b;
1060e1cf1444SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1061e1cf1444SMark Adams }
1062e1cf1444SMark Adams 
1063e1cf1444SMark Adams /*@
1064e1cf1444SMark Adams   PCGAMGSetRecomputeEstEig - Set flag for Chebyshev smoothers to recompute the eigen estimates
1065e1cf1444SMark Adams 
1066e1cf1444SMark Adams   Collective
1067e1cf1444SMark Adams 
1068e1cf1444SMark Adams   Input Parameters:
1069e1cf1444SMark Adams + pc - the preconditioner context
1070e1cf1444SMark Adams - b  - flag
1071e1cf1444SMark Adams 
1072e1cf1444SMark Adams   Options Database Key:
1073e1cf1444SMark Adams . -pc_gamg_recompute_esteig <true> - use the eigen estimate
1074e1cf1444SMark Adams 
1075e1cf1444SMark Adams   Level: advanced
1076e1cf1444SMark Adams 
1077562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
1078e1cf1444SMark Adams @*/
1079e1cf1444SMark Adams PetscErrorCode PCGAMGSetRecomputeEstEig(PC pc, PetscBool b)
1080e1cf1444SMark Adams {
1081e1cf1444SMark Adams   PetscFunctionBegin;
1082e1cf1444SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1083e1cf1444SMark Adams   PetscTryMethod(pc, "PCGAMGSetRecomputeEstEig_C", (PC, PetscBool), (pc, b));
1084e1cf1444SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1085e1cf1444SMark Adams }
1086e1cf1444SMark Adams 
1087e1cf1444SMark Adams static PetscErrorCode PCGAMGSetRecomputeEstEig_GAMG(PC pc, PetscBool b)
1088e1cf1444SMark Adams {
1089e1cf1444SMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1090e1cf1444SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1091e1cf1444SMark Adams 
1092e1cf1444SMark Adams   PetscFunctionBegin;
1093e1cf1444SMark Adams   pc_gamg->recompute_esteig = b;
10943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109518c3aa7eSMark }
109618c3aa7eSMark 
109718c3aa7eSMark /*@
1098f1580f4eSBarry Smith   PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY?
109918c3aa7eSMark 
1100c3339decSBarry Smith   Collective
110118c3aa7eSMark 
110218c3aa7eSMark   Input Parameters:
110318c3aa7eSMark + pc   - the preconditioner context
1104feefa0e1SJacob Faibussowitsch . emax - max eigenvalue
1105feefa0e1SJacob Faibussowitsch - emin - min eigenvalue
110618c3aa7eSMark 
110718c3aa7eSMark   Options Database Key:
1108147403d9SBarry Smith . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
110918c3aa7eSMark 
111018c3aa7eSMark   Level: intermediate
111118c3aa7eSMark 
1112562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetUseSAEstEig()`
111318c3aa7eSMark @*/
1114d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin)
1115d71ae5a4SJacob Faibussowitsch {
111618c3aa7eSMark   PetscFunctionBegin;
111718c3aa7eSMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1118cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin));
11193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112018c3aa7eSMark }
112141ffd417SStefano Zampini 
1122d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin)
1123d71ae5a4SJacob Faibussowitsch {
112418c3aa7eSMark   PC_MG   *mg      = (PC_MG *)pc->data;
112518c3aa7eSMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
112618c3aa7eSMark 
112718c3aa7eSMark   PetscFunctionBegin;
11285f80ce2aSJacob 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);
11295f80ce2aSJacob 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);
113018c3aa7eSMark   pc_gamg->emax = emax;
113118c3aa7eSMark   pc_gamg->emin = emin;
11323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113318c3aa7eSMark }
113418c3aa7eSMark 
113518c3aa7eSMark /*@
1136f1580f4eSBarry Smith   PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner
1137dfd5c07aSMark F. Adams 
1138c3339decSBarry Smith   Collective
1139dfd5c07aSMark F. Adams 
1140dfd5c07aSMark F. Adams   Input Parameters:
11411cc46a46SBarry Smith + pc - the preconditioner context
1142f1580f4eSBarry Smith - n  - `PETSC_TRUE` or `PETSC_FALSE`
1143dfd5c07aSMark F. Adams 
1144dfd5c07aSMark F. Adams   Options Database Key:
1145147403d9SBarry Smith . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1146dfd5c07aSMark F. Adams 
1147dfd5c07aSMark F. Adams   Level: intermediate
1148dfd5c07aSMark F. Adams 
1149f1580f4eSBarry Smith   Note:
1150147403d9SBarry Smith   May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1151cab9ed1eSBarry Smith   rebuilding the preconditioner quicker.
1152cab9ed1eSBarry Smith 
1153562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`
1154dfd5c07aSMark F. Adams @*/
1155d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1156d71ae5a4SJacob Faibussowitsch {
1157dfd5c07aSMark F. Adams   PetscFunctionBegin;
1158dfd5c07aSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1159cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n));
11603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1161dfd5c07aSMark F. Adams }
1162dfd5c07aSMark F. Adams 
1163d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1164d71ae5a4SJacob Faibussowitsch {
1165dfd5c07aSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1166dfd5c07aSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1167dfd5c07aSMark F. Adams 
1168dfd5c07aSMark F. Adams   PetscFunctionBegin;
1169dfd5c07aSMark F. Adams   pc_gamg->reuse_prol = n;
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1171dfd5c07aSMark F. Adams }
1172dfd5c07aSMark F. Adams 
1173ffc955d6SMark F. Adams /*@
1174f1580f4eSBarry 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
1175f1580f4eSBarry Smith   used as the smoother
1176ffc955d6SMark F. Adams 
1177c3339decSBarry Smith   Collective
1178ffc955d6SMark F. Adams 
1179ffc955d6SMark F. Adams   Input Parameters:
1180cab9ed1eSBarry Smith + pc  - the preconditioner context
1181f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not
1182ffc955d6SMark F. Adams 
1183ffc955d6SMark F. Adams   Options Database Key:
1184147403d9SBarry Smith . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1185ffc955d6SMark F. Adams 
1186ffc955d6SMark F. Adams   Level: intermediate
1187ffc955d6SMark F. Adams 
1188562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCASM`, `PCSetType`
1189ffc955d6SMark F. Adams @*/
1190d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1191d71ae5a4SJacob Faibussowitsch {
1192ffc955d6SMark F. Adams   PetscFunctionBegin;
1193ffc955d6SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1194cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg));
11953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1196ffc955d6SMark F. Adams }
1197ffc955d6SMark F. Adams 
1198d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1199d71ae5a4SJacob Faibussowitsch {
1200ffc955d6SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1201ffc955d6SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1202ffc955d6SMark F. Adams 
1203ffc955d6SMark F. Adams   PetscFunctionBegin;
1204cab9ed1eSBarry Smith   pc_gamg->use_aggs_in_asm = flg;
12053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1206ffc955d6SMark F. Adams }
1207ffc955d6SMark F. Adams 
1208171cca9aSMark Adams /*@
1209d529f056Smarkadams4   PCGAMGSetParallelCoarseGridSolve - allow a parallel coarse grid solver
1210171cca9aSMark Adams 
1211c3339decSBarry Smith   Collective
1212171cca9aSMark Adams 
1213171cca9aSMark Adams   Input Parameters:
1214171cca9aSMark Adams + pc  - the preconditioner context
1215f1580f4eSBarry Smith - flg - `PETSC_TRUE` to not force coarse grid onto one processor
1216171cca9aSMark Adams 
1217171cca9aSMark Adams   Options Database Key:
1218d529f056Smarkadams4 . -pc_gamg_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1219171cca9aSMark Adams 
1220171cca9aSMark Adams   Level: intermediate
1221171cca9aSMark Adams 
1222562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`
1223171cca9aSMark Adams @*/
1224d529f056Smarkadams4 PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC pc, PetscBool flg)
1225d71ae5a4SJacob Faibussowitsch {
1226171cca9aSMark Adams   PetscFunctionBegin;
1227171cca9aSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1228d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGSetParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg));
12293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1230171cca9aSMark Adams }
1231171cca9aSMark Adams 
1232d529f056Smarkadams4 static PetscErrorCode PCGAMGSetParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1233d71ae5a4SJacob Faibussowitsch {
1234171cca9aSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1235171cca9aSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1236171cca9aSMark Adams 
1237171cca9aSMark Adams   PetscFunctionBegin;
1238171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = flg;
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1240ffc955d6SMark F. Adams }
1241ffc955d6SMark F. Adams 
12424ef23d27SMark F. Adams /*@
1243f1580f4eSBarry 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
1244ce7c7f2fSMark Adams 
1245c3339decSBarry Smith   Collective
1246ce7c7f2fSMark Adams 
1247ce7c7f2fSMark Adams   Input Parameters:
1248ce7c7f2fSMark Adams + pc  - the preconditioner context
1249f1580f4eSBarry Smith - flg - `PETSC_TRUE` to pin coarse grids to the CPU
1250ce7c7f2fSMark Adams 
1251ce7c7f2fSMark Adams   Options Database Key:
1252147403d9SBarry Smith . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1253ce7c7f2fSMark Adams 
1254ce7c7f2fSMark Adams   Level: intermediate
1255ce7c7f2fSMark Adams 
1256562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetParallelCoarseGridSolve()`
1257ce7c7f2fSMark Adams @*/
1258d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1259d71ae5a4SJacob Faibussowitsch {
1260ce7c7f2fSMark Adams   PetscFunctionBegin;
1261ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1262cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg));
12633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1264ce7c7f2fSMark Adams }
1265ce7c7f2fSMark Adams 
1266d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1267d71ae5a4SJacob Faibussowitsch {
1268ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1269ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1270ce7c7f2fSMark Adams 
1271ce7c7f2fSMark Adams   PetscFunctionBegin;
1272ce7c7f2fSMark Adams   pc_gamg->cpu_pin_coarse_grids = flg;
12733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1274ce7c7f2fSMark Adams }
1275ce7c7f2fSMark Adams 
1276ce7c7f2fSMark Adams /*@
1277147403d9SBarry Smith   PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1278ce7c7f2fSMark Adams 
1279c3339decSBarry Smith   Collective
1280ce7c7f2fSMark Adams 
1281ce7c7f2fSMark Adams   Input Parameters:
1282ce7c7f2fSMark Adams + pc  - the preconditioner context
1283f1580f4eSBarry Smith - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD`
1284ce7c7f2fSMark Adams 
1285ce7c7f2fSMark Adams   Options Database Key:
1286147403d9SBarry Smith . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1287ce7c7f2fSMark Adams 
1288ce7c7f2fSMark Adams   Level: intermediate
1289ce7c7f2fSMark Adams 
1290562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD`
1291ce7c7f2fSMark Adams @*/
1292d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1293d71ae5a4SJacob Faibussowitsch {
1294ce7c7f2fSMark Adams   PetscFunctionBegin;
1295ce7c7f2fSMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1296cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg));
12973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1298ce7c7f2fSMark Adams }
1299ce7c7f2fSMark Adams 
1300d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1301d71ae5a4SJacob Faibussowitsch {
1302ce7c7f2fSMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1303ce7c7f2fSMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1304ce7c7f2fSMark Adams 
1305ce7c7f2fSMark Adams   PetscFunctionBegin;
1306ce7c7f2fSMark Adams   pc_gamg->layout_type = flg;
13073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1308ce7c7f2fSMark Adams }
1309ce7c7f2fSMark Adams 
1310ce7c7f2fSMark Adams /*@
1311f1580f4eSBarry Smith   PCGAMGSetNlevels -  Sets the maximum number of levels `PCGAMG` will use
13124ef23d27SMark F. Adams 
1313*8926f930SMark Adams   Collective
13144ef23d27SMark F. Adams 
13154ef23d27SMark F. Adams   Input Parameters:
13161cc46a46SBarry Smith + pc - the preconditioner
13171cc46a46SBarry Smith - n  - the maximum number of levels to use
13184ef23d27SMark F. Adams 
13194ef23d27SMark F. Adams   Options Database Key:
1320147403d9SBarry Smith . -pc_mg_levels <n> - set the maximum number of levels to allow
13214ef23d27SMark F. Adams 
13224ef23d27SMark F. Adams   Level: intermediate
13234ef23d27SMark F. Adams 
1324feefa0e1SJacob Faibussowitsch   Developer Notes:
1325f1580f4eSBarry Smith   Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`
1326f1580f4eSBarry Smith 
1327562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`
13284ef23d27SMark F. Adams @*/
1329d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1330d71ae5a4SJacob Faibussowitsch {
13314ef23d27SMark F. Adams   PetscFunctionBegin;
13324ef23d27SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1333cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n));
13343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13354ef23d27SMark F. Adams }
13364ef23d27SMark F. Adams 
1337d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1338d71ae5a4SJacob Faibussowitsch {
13394ef23d27SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
13404ef23d27SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
13414ef23d27SMark F. Adams 
13424ef23d27SMark F. Adams   PetscFunctionBegin;
13439d5b6da9SMark F. Adams   pc_gamg->Nlevels = n;
13443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13454ef23d27SMark F. Adams }
13464ef23d27SMark F. Adams 
13473542efc5SMark F. Adams /*@
1348*8926f930SMark Adams   PCGAMGASMSetHEM -  Sets the number of HEM matching passed
1349*8926f930SMark Adams 
1350*8926f930SMark Adams   Collective
1351*8926f930SMark Adams 
1352*8926f930SMark Adams   Input Parameters:
1353*8926f930SMark Adams + pc - the preconditioner
1354*8926f930SMark Adams - n  - number of HEM matching passed to construct ASM subdomains
1355*8926f930SMark Adams 
1356*8926f930SMark Adams   Options Database Key:
1357*8926f930SMark Adams . -pc_gamg_asm_hem <n> - set the number of HEM matching passed
1358*8926f930SMark Adams 
1359*8926f930SMark Adams   Level: intermediate
1360*8926f930SMark Adams 
1361*8926f930SMark Adams   Developer Notes:
1362*8926f930SMark Adams   Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`
1363*8926f930SMark Adams 
1364*8926f930SMark Adams .seealso: [](ch_ksp), `PCGAMG`
1365*8926f930SMark Adams @*/
1366*8926f930SMark Adams PetscErrorCode PCGAMGASMSetHEM(PC pc, PetscInt n)
1367*8926f930SMark Adams {
1368*8926f930SMark Adams   PetscFunctionBegin;
1369*8926f930SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1370*8926f930SMark Adams   PetscTryMethod(pc, "PCGAMGASMSetHEM_C", (PC, PetscInt), (pc, n));
1371*8926f930SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1372*8926f930SMark Adams }
1373*8926f930SMark Adams 
1374*8926f930SMark Adams static PetscErrorCode PCGAMGASMSetHEM_GAMG(PC pc, PetscInt n)
1375*8926f930SMark Adams {
1376*8926f930SMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1377*8926f930SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1378*8926f930SMark Adams 
1379*8926f930SMark Adams   PetscFunctionBegin;
1380*8926f930SMark Adams   pc_gamg->asm_hem_aggs = n;
1381*8926f930SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1382*8926f930SMark Adams }
1383*8926f930SMark Adams 
1384*8926f930SMark Adams /*@
13853542efc5SMark F. Adams   PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
13863542efc5SMark F. Adams 
138720f4b53cSBarry Smith   Not Collective
13883542efc5SMark F. Adams 
13893542efc5SMark F. Adams   Input Parameters:
13901cc46a46SBarry Smith + pc - the preconditioner context
1391c9567895SMark . 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
1392055c8bd0SJed Brown - n  - number of threshold values provided in array
13933542efc5SMark F. Adams 
13943542efc5SMark F. Adams   Options Database Key:
1395147403d9SBarry Smith . -pc_gamg_threshold <threshold> - the threshold to drop edges
13963542efc5SMark F. Adams 
139720f4b53cSBarry Smith   Level: intermediate
139820f4b53cSBarry Smith 
139995452b02SPatrick Sanan   Notes:
1400af3c827dSMark 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.
1401f1580f4eSBarry 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.
1402cab9ed1eSBarry Smith 
140320f4b53cSBarry Smith   If `n` is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening.
1404f1580f4eSBarry Smith   In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`.
140520f4b53cSBarry Smith   If `n` is greater than the total number of levels, the excess entries in threshold will not be used.
14063542efc5SMark F. Adams 
1407562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetMinDegreeOrderingMISk()`, `PCGAMGSetThresholdScale()`
14083542efc5SMark F. Adams @*/
1409d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1410d71ae5a4SJacob Faibussowitsch {
14113542efc5SMark F. Adams   PetscFunctionBegin;
14123542efc5SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14134f572ea9SToby Isaac   if (n) PetscAssertPointer(v, 2);
1414cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n));
14153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14163542efc5SMark F. Adams }
14173542efc5SMark F. Adams 
1418d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1419d71ae5a4SJacob Faibussowitsch {
1420c20e4228SMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1421c20e4228SMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1422c1eae691SMark Adams   PetscInt i;
1423c1eae691SMark Adams   PetscFunctionBegin;
1424055c8bd0SJed Brown   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1425055c8bd0SJed Brown   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
14263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1427c1eae691SMark Adams }
1428c1eae691SMark Adams 
1429c1eae691SMark Adams /*@
1430f1580f4eSBarry Smith   PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids
1431c9567895SMark 
1432c3339decSBarry Smith   Collective
1433c9567895SMark 
1434c9567895SMark   Input Parameters:
1435c9567895SMark + pc - the preconditioner context
1436f1580f4eSBarry Smith . v  - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA
1437c9567895SMark - n  - number of values provided in array
1438c9567895SMark 
1439c9567895SMark   Options Database Key:
1440147403d9SBarry Smith . -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1441c9567895SMark 
1442c9567895SMark   Level: intermediate
1443c9567895SMark 
1444562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`
1445c9567895SMark @*/
1446d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1447d71ae5a4SJacob Faibussowitsch {
1448c9567895SMark   PetscFunctionBegin;
1449c9567895SMark   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14504f572ea9SToby Isaac   if (n) PetscAssertPointer(v, 2);
1451cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n));
14523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1453c9567895SMark }
1454c9567895SMark 
1455d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1456d71ae5a4SJacob Faibussowitsch {
1457c9567895SMark   PC_MG   *mg      = (PC_MG *)pc->data;
1458c9567895SMark   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1459c9567895SMark   PetscInt i;
1460c9567895SMark   PetscFunctionBegin;
1461c9567895SMark   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1462c9567895SMark   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
14633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1464c9567895SMark }
1465c9567895SMark 
1466c9567895SMark /*@
1467c1eae691SMark Adams   PCGAMGSetThresholdScale - Relative threshold reduction at each level
1468c1eae691SMark Adams 
146920f4b53cSBarry Smith   Not Collective
1470c1eae691SMark Adams 
1471c1eae691SMark Adams   Input Parameters:
1472c1eae691SMark Adams + pc - the preconditioner context
1473feefa0e1SJacob Faibussowitsch - v  - the threshold value reduction, usually < 1.0
1474c1eae691SMark Adams 
1475c1eae691SMark Adams   Options Database Key:
1476147403d9SBarry Smith . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1477c1eae691SMark Adams 
147820f4b53cSBarry Smith   Level: advanced
147920f4b53cSBarry Smith 
1480f1580f4eSBarry Smith   Note:
1481f1580f4eSBarry Smith   The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`.
1482f1580f4eSBarry Smith   This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`.
1483055c8bd0SJed Brown 
1484562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`
1485c1eae691SMark Adams @*/
1486d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1487d71ae5a4SJacob Faibussowitsch {
14883542efc5SMark F. Adams   PetscFunctionBegin;
1489c1eae691SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1490cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v));
14913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1492c1eae691SMark Adams }
1493c1eae691SMark Adams 
1494d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1495d71ae5a4SJacob Faibussowitsch {
1496c1eae691SMark Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1497c1eae691SMark Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1498c1eae691SMark Adams   PetscFunctionBegin;
1499c1eae691SMark Adams   pc_gamg->threshold_scale = v;
15003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15013542efc5SMark F. Adams }
15023542efc5SMark F. Adams 
1503e20c40e8SBarry Smith /*@C
1504f1580f4eSBarry Smith   PCGAMGSetType - Set the type of algorithm `PCGAMG` should use
1505676e1743SMark F. Adams 
1506c3339decSBarry Smith   Collective
1507676e1743SMark F. Adams 
1508676e1743SMark F. Adams   Input Parameters:
1509c60c7ad4SBarry Smith + pc   - the preconditioner context
1510f1580f4eSBarry Smith - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL`
1511676e1743SMark F. Adams 
1512676e1743SMark F. Adams   Options Database Key:
15135e4ac8c8Smarkadams4 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply - only agg supported
1514676e1743SMark F. Adams 
1515676e1743SMark F. Adams   Level: intermediate
1516676e1743SMark F. Adams 
1517562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType`
1518676e1743SMark F. Adams @*/
1519d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1520d71ae5a4SJacob Faibussowitsch {
1521676e1743SMark F. Adams   PetscFunctionBegin;
1522676e1743SMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1523cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type));
15243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1525676e1743SMark F. Adams }
1526676e1743SMark F. Adams 
1527e20c40e8SBarry Smith /*@C
1528f1580f4eSBarry Smith   PCGAMGGetType - Get the type of algorithm `PCGAMG` will use
1529c60c7ad4SBarry Smith 
1530c3339decSBarry Smith   Collective
1531c60c7ad4SBarry Smith 
1532c60c7ad4SBarry Smith   Input Parameter:
1533c60c7ad4SBarry Smith . pc - the preconditioner context
1534c60c7ad4SBarry Smith 
1535c60c7ad4SBarry Smith   Output Parameter:
1536c60c7ad4SBarry Smith . type - the type of algorithm used
1537c60c7ad4SBarry Smith 
1538c60c7ad4SBarry Smith   Level: intermediate
1539c60c7ad4SBarry Smith 
1540562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType`
1541c60c7ad4SBarry Smith @*/
1542d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1543d71ae5a4SJacob Faibussowitsch {
1544c60c7ad4SBarry Smith   PetscFunctionBegin;
1545c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1546cac4c232SBarry Smith   PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type));
15473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1548c60c7ad4SBarry Smith }
1549c60c7ad4SBarry Smith 
1550d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1551d71ae5a4SJacob Faibussowitsch {
1552c60c7ad4SBarry Smith   PC_MG   *mg      = (PC_MG *)pc->data;
1553c60c7ad4SBarry Smith   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1554c60c7ad4SBarry Smith 
1555c60c7ad4SBarry Smith   PetscFunctionBegin;
1556c60c7ad4SBarry Smith   *type = pc_gamg->type;
15573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1558c60c7ad4SBarry Smith }
1559c60c7ad4SBarry Smith 
1560d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1561d71ae5a4SJacob Faibussowitsch {
15621ab5ffc9SJed Brown   PC_MG   *mg      = (PC_MG *)pc->data;
15631ab5ffc9SJed Brown   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
15645f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PC);
1565676e1743SMark F. Adams 
1566676e1743SMark F. Adams   PetscFunctionBegin;
1567c60c7ad4SBarry Smith   pc_gamg->type = type;
15689566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(GAMGList, type, &r));
15696adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type);
15701ab5ffc9SJed Brown   if (pc_gamg->ops->destroy) {
15719566063dSJacob Faibussowitsch     PetscCall((*pc_gamg->ops->destroy)(pc));
15729566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps)));
1573e616c208SToby Isaac     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1574da81f932SPierre Jolivet     /* cleaning up common data in pc_gamg - this should disappear someday */
15753ae0bb68SMark Adams     pc_gamg->data_cell_cols      = 0;
15763ae0bb68SMark Adams     pc_gamg->data_cell_rows      = 0;
15773ae0bb68SMark Adams     pc_gamg->orig_data_cell_cols = 0;
15783ae0bb68SMark Adams     pc_gamg->orig_data_cell_rows = 0;
15799566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
15803ae0bb68SMark Adams     pc_gamg->data_sz = 0;
15811ab5ffc9SJed Brown   }
15829566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->gamg_type_name));
15839566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name));
15849566063dSJacob Faibussowitsch   PetscCall((*r)(pc));
15853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1586676e1743SMark F. Adams }
1587676e1743SMark F. Adams 
1588d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer)
1589d71ae5a4SJacob Faibussowitsch {
15905adeb434SBarry Smith   PC_MG    *mg      = (PC_MG *)pc->data;
15915adeb434SBarry Smith   PC_GAMG  *pc_gamg = (PC_GAMG *)mg->innerctx;
1592e7d4b4cbSMark Adams   PetscReal gc = 0, oc = 0;
159390db8557SMark Adams 
15945adeb434SBarry Smith   PetscFunctionBegin;
15959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "    GAMG specific options\n"));
15969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold for dropping small values in graph on each level ="));
15979566063dSJacob Faibussowitsch   for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i]));
15989566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
15999566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale));
1600*8926f930SMark Adams   if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates from coarsening process to define subdomains for PCASM\n")); // this take presedence
1601*8926f930SMark Adams   else if (pc_gamg->asm_hem_aggs) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates made with %d applications of heavy edge matching (HEM) to define subdomains for PCASM\n", (int)pc_gamg->asm_hem_aggs));
160248a46eb9SPierre 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"));
16031baa6e33SBarry Smith   if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer));
16049566063dSJacob Faibussowitsch   PetscCall(PCMGGetGridComplexity(pc, &gc, &oc));
160563a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      Complexity:    grid = %g    operator = %g\n", (double)gc, (double)oc));
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16075adeb434SBarry Smith }
16085adeb434SBarry Smith 
160966976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject)
1610d71ae5a4SJacob Faibussowitsch {
1611676e1743SMark F. Adams   PC_MG             *mg      = (PC_MG *)pc->data;
1612676e1743SMark F. Adams   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
16137e6512fdSJed Brown   PetscBool          flag;
16143b4367a7SBarry Smith   MPI_Comm           comm;
161518c3aa7eSMark   char               prefix[256], tname[32];
1616c1eae691SMark Adams   PetscInt           i, n;
161714a9496bSBarry Smith   const char        *pcpre;
16180a545947SLisandro Dalcin   static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL};
1619*8926f930SMark Adams 
16205b89ad90SMark F. Adams   PetscFunctionBegin;
16219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1622d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options");
16235e4ac8c8Smarkadams4   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));
16241baa6e33SBarry Smith   if (flag) PetscCall(PCGAMGSetType(pc, tname));
16259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL));
1626f1580f4eSBarry 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));
1627e1cf1444SMark 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));
16289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL));
16299566063dSJacob 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));
1630d529f056Smarkadams4   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));
16319566063dSJacob 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));
16329371c9d4SSatish 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,
16339371c9d4SSatish Balay                              (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL));
16349566063dSJacob 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));
16359566063dSJacob 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));
1636*8926f930SMark Adams   PetscCall(PetscOptionsInt("-pc_gamg_asm_hem_aggs", "Number of HEM matching passed in aggregates for ASM smoother", "PCGAMGASMSetHEM", pc_gamg->asm_hem_aggs, &pc_gamg->asm_hem_aggs, NULL));
16379566063dSJacob 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));
163818c3aa7eSMark   n = PETSC_MG_MAXLEVELS;
16399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag));
164018c3aa7eSMark   if (!flag || n < PETSC_MG_MAXLEVELS) {
1641efd3c5ceSMark Adams     if (!flag) n = 1;
1642c1eae691SMark Adams     i = n;
1643d71ae5a4SJacob Faibussowitsch     do {
1644d71ae5a4SJacob Faibussowitsch       pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1645d71ae5a4SJacob Faibussowitsch     } while (++i < PETSC_MG_MAXLEVELS);
1646c1eae691SMark Adams   }
16479c15c1aeSMark Adams   PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels (should get from base class)", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL));
16489c15c1aeSMark Adams   PetscCheck(pc_gamg->Nlevels <= PETSC_MG_MAXLEVELS, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_mg_levels (%d) >= PETSC_MG_MAXLEVELS (%d)", (int)pc_gamg->Nlevels, (int)PETSC_MG_MAXLEVELS);
1649c9567895SMark   n = PETSC_MG_MAXLEVELS;
16509566063dSJacob 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));
1651c9567895SMark   if (!flag) i = 0;
1652c9567895SMark   else i = n;
1653d71ae5a4SJacob Faibussowitsch   do {
1654d71ae5a4SJacob Faibussowitsch     pc_gamg->level_reduction_factors[i] = -1;
1655d71ae5a4SJacob Faibussowitsch   } while (++i < PETSC_MG_MAXLEVELS);
165618c3aa7eSMark   {
165718c3aa7eSMark     PetscReal eminmax[2] = {0., 0.};
165818c3aa7eSMark     n                    = 2;
16599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag));
166018c3aa7eSMark     if (flag) {
166108401ef6SPierre Jolivet       PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
16629566063dSJacob Faibussowitsch       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
166318c3aa7eSMark     }
166418c3aa7eSMark   }
16659c15c1aeSMark Adams 
1666b7cbab4eSMark Adams   /* set options for subtype */
1667dbbe0bcdSBarry Smith   PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject));
166818c3aa7eSMark 
16699566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
16709566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
1671d0609cedSBarry Smith   PetscOptionsHeadEnd();
16723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16735b89ad90SMark F. Adams }
16745b89ad90SMark F. Adams 
16755b89ad90SMark F. Adams /*MC
16761cc46a46SBarry Smith   PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
16775b89ad90SMark F. Adams 
1678280d9858SJed Brown   Options Database Keys:
16795e4ac8c8Smarkadams4 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical (only smoothed aggregation, agg, supported)
1680da81f932SPierre Jolivet . -pc_gamg_repartition  <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined
168121d928e4Smarkadams4 . -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
1682da81f932SPierre 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>
1683cab9ed1eSBarry Smith                                         equations on each process that has degrees of freedom
16842d776b49SBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
168521d928e4Smarkadams4 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true)
168621d928e4Smarkadams4 . -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)
16872d776b49SBarry Smith - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1688cab9ed1eSBarry Smith 
1689f1580f4eSBarry Smith   Options Database Keys for Aggregation:
1690cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1691d529f056Smarkadams4 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1692d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1693d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1694*8926f930SMark Adams . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1695d529f056Smarkadams4 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1696cab9ed1eSBarry Smith 
1697f1580f4eSBarry Smith   Options Database Keys for Multigrid:
1698a9f5add0SYANG Zongze + -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()`
1699db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1700db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
170121d928e4Smarkadams4 - -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
17025b89ad90SMark F. Adams 
170320f4b53cSBarry Smith   Level: intermediate
170420f4b53cSBarry Smith 
170595452b02SPatrick Sanan   Notes:
1706f1580f4eSBarry Smith   To obtain good performance for `PCGAMG` for vector valued problems you must
1707f1580f4eSBarry Smith   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point
1708f1580f4eSBarry Smith   call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1709f1580f4eSBarry Smith 
171004c3f3b8SBarry Smith   The many options for `PCMG` also work directly for `PCGAMG` such as controlling the smoothers on each level etc.
171104c3f3b8SBarry Smith 
171204c3f3b8SBarry Smith   See [the Users Manual section on PCGAMG](sec_amg) and [the Users Manual section on PCMG](sec_mg)for more details.
17131cc46a46SBarry Smith 
1714562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`,
1715d529f056Smarkadams4           `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
17165b89ad90SMark F. Adams M*/
1717d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1718d71ae5a4SJacob Faibussowitsch {
17195b89ad90SMark F. Adams   PC_GAMG *pc_gamg;
17205b89ad90SMark F. Adams   PC_MG   *mg;
17215b89ad90SMark F. Adams 
17225b89ad90SMark F. Adams   PetscFunctionBegin;
17231c1aac46SBarry Smith   /* register AMG type */
17249566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
17251c1aac46SBarry Smith 
17265b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
17279566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
17289566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
17295b89ad90SMark F. Adams 
17305b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
17314dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg));
17329566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
17335b89ad90SMark F. Adams   mg           = (PC_MG *)pc->data;
17345b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
17355b89ad90SMark F. Adams 
17364dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg->ops));
17371ab5ffc9SJed Brown 
17389d5b6da9SMark F. Adams   /* these should be in subctx but repartitioning needs simple arrays */
17399d5b6da9SMark F. Adams   pc_gamg->data_sz = 0;
17400a545947SLisandro Dalcin   pc_gamg->data    = NULL;
17415b89ad90SMark F. Adams 
17429d5b6da9SMark F. Adams   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
17435b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
17445b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
17455b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
17465b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
17475adeb434SBarry Smith   mg->view                = PCView_GAMG;
17485b89ad90SMark F. Adams 
17499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG));
17509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG));
17519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG));
17529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG));
17539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG));
1754e1cf1444SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", PCGAMGSetRecomputeEstEig_GAMG));
17559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG));
17569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG));
1757d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", PCGAMGSetParallelCoarseGridSolve_GAMG));
17589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG));
17599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG));
17609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG));
17619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG));
17629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG));
17639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG));
17649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG));
17659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG));
1766*8926f930SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", PCGAMGASMSetHEM_GAMG));
17679d5b6da9SMark F. Adams   pc_gamg->repart                          = PETSC_FALSE;
176821d928e4Smarkadams4   pc_gamg->reuse_prol                      = PETSC_TRUE;
17690c3bc534SBarry Smith   pc_gamg->use_aggs_in_asm                 = PETSC_FALSE;
1770171cca9aSMark Adams   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1771a0095786SMark   pc_gamg->cpu_pin_coarse_grids            = PETSC_FALSE;
1772a0095786SMark   pc_gamg->layout_type                     = PCGAMG_LAYOUT_SPREAD;
1773038f3aa4SMark F. Adams   pc_gamg->min_eq_proc                     = 50;
1774*8926f930SMark Adams   pc_gamg->asm_hem_aggs                    = 0;
177525a145a7SMark Adams   pc_gamg->coarse_eq_limit                 = 50;
177653134ebeSMark Adams   for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1;
1777c1eae691SMark Adams   pc_gamg->threshold_scale  = 1.;
177818c3aa7eSMark   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
17799ab59c8bSMark Adams   pc_gamg->current_level    = 0; /* don't need to init really */
17807e6512fdSJed Brown   pc_gamg->use_sa_esteig    = PETSC_TRUE;
1781e1cf1444SMark Adams   pc_gamg->recompute_esteig = PETSC_TRUE;
178218c3aa7eSMark   pc_gamg->emin             = 0;
178318c3aa7eSMark   pc_gamg->emax             = 0;
178418c3aa7eSMark 
1785c238b0ebSToby Isaac   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
17869d5b6da9SMark F. Adams 
1787bd94a7aaSJed Brown   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
17889566063dSJacob Faibussowitsch   PetscCall(PCGAMGSetType(pc, PCGAMGAGG));
17893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17905b89ad90SMark F. Adams }
17913e3471ccSMark Adams 
17923e3471ccSMark Adams /*@C
1793f1580f4eSBarry Smith   PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called
1794f1580f4eSBarry Smith   from `PCInitializePackage()`.
17953e3471ccSMark Adams 
17963e3471ccSMark Adams   Level: developer
17973e3471ccSMark Adams 
1798562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscInitialize()`
17993e3471ccSMark Adams @*/
1800d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGInitializePackage(void)
1801d71ae5a4SJacob Faibussowitsch {
18024555aa8cSStefano Zampini   PetscInt l;
18033e3471ccSMark Adams 
18043e3471ccSMark Adams   PetscFunctionBegin;
18053ba16761SJacob Faibussowitsch   if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
18063e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_TRUE;
18079566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO));
18089566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG));
18099566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical));
18109566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1811c1c463dbSMark Adams 
1812c1c463dbSMark Adams   /* general events */
1813849bee69SMark Adams   PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]));
1814849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]));
1815849bee69SMark Adams   PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]));
1816849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]));
1817849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]));
1818849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]));
1819849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]));
1820849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]));
1821849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]));
1822849bee69SMark Adams   PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]));
1823849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]));
1824849bee69SMark Adams   PetscCall(PetscLogEventRegister("  GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]));
1825849bee69SMark Adams   PetscCall(PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]));
1826849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */
1827849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */
1828849bee69SMark Adams   /* PetscCall(PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */
18294555aa8cSStefano Zampini   for (l = 0; l < PETSC_MG_MAXLEVELS; l++) {
18304555aa8cSStefano Zampini     char ename[32];
18315b89ad90SMark F. Adams 
183263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l));
18339566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
183463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l));
18359566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
183663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l));
18379566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
18384555aa8cSStefano Zampini   }
18394555aa8cSStefano Zampini #if defined(GAMG_STAGES)
1840849bee69SMark Adams   { /* create timer stages */
18415b89ad90SMark F. Adams     char str[32];
1842a364092eSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0));
18439566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
18445b89ad90SMark F. Adams   }
18455b89ad90SMark F. Adams #endif
18463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18473e3471ccSMark Adams }
18483e3471ccSMark Adams 
18493e3471ccSMark Adams /*@C
1850f1580f4eSBarry Smith   PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is
1851f1580f4eSBarry Smith   called from `PetscFinalize()` automatically.
18523e3471ccSMark Adams 
18533e3471ccSMark Adams   Level: developer
18543e3471ccSMark Adams 
1855562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscFinalize()`
18563e3471ccSMark Adams @*/
1857d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGFinalizePackage(void)
1858d71ae5a4SJacob Faibussowitsch {
18593e3471ccSMark Adams   PetscFunctionBegin;
18603e3471ccSMark Adams   PCGAMGPackageInitialized = PETSC_FALSE;
18619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&GAMGList));
18623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18633e3471ccSMark Adams }
1864a36cf38bSToby Isaac 
1865a36cf38bSToby Isaac /*@C
1866f1580f4eSBarry Smith   PCGAMGRegister - Register a `PCGAMG` implementation.
1867a36cf38bSToby Isaac 
1868a36cf38bSToby Isaac   Input Parameters:
1869f1580f4eSBarry Smith + type   - string that will be used as the name of the `PCGAMG` type.
1870a36cf38bSToby Isaac - create - function for creating the gamg context.
1871a36cf38bSToby Isaac 
1872f1580f4eSBarry Smith   Level: developer
1873a36cf38bSToby Isaac 
1874562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1875a36cf38bSToby Isaac @*/
1876d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1877d71ae5a4SJacob Faibussowitsch {
1878a36cf38bSToby Isaac   PetscFunctionBegin;
18799566063dSJacob Faibussowitsch   PetscCall(PCGAMGInitializePackage());
18809566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&GAMGList, type, create));
18813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1882a36cf38bSToby Isaac }
18832d776b49SBarry Smith 
18842d776b49SBarry Smith /*@C
18852d776b49SBarry Smith   PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process
18862d776b49SBarry Smith 
18872d776b49SBarry Smith   Input Parameters:
18882d776b49SBarry Smith + pc - the `PCGAMG`
18892d776b49SBarry Smith - A  - the matrix, for any level
18902d776b49SBarry Smith 
18912d776b49SBarry Smith   Output Parameter:
18922d776b49SBarry Smith . G - the graph
18932d776b49SBarry Smith 
18942d776b49SBarry Smith   Level: advanced
18952d776b49SBarry Smith 
1896562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
18972d776b49SBarry Smith @*/
18982d776b49SBarry Smith PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G)
18992d776b49SBarry Smith {
19002d776b49SBarry Smith   PC_MG   *mg      = (PC_MG *)pc->data;
19012d776b49SBarry Smith   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
19022d776b49SBarry Smith 
19032d776b49SBarry Smith   PetscFunctionBegin;
19042d776b49SBarry Smith   PetscCall(pc_gamg->ops->creategraph(pc, A, G));
19053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19062d776b49SBarry Smith }
1907