1e6f8f311SMark Adams /*
2e6f8f311SMark Adams This file implements an AmgX preconditioner in PETSc as part of PC.
3e6f8f311SMark Adams */
4e6f8f311SMark Adams
5e6f8f311SMark Adams /*
6e6f8f311SMark Adams Include files needed for the AmgX preconditioner:
7e6f8f311SMark Adams pcimpl.h - private include file intended for use by all preconditioners
8e6f8f311SMark Adams */
9e6f8f311SMark Adams
10e6f8f311SMark Adams #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
110e6b6b59SJacob Faibussowitsch #include <petscdevice_cuda.h>
12e6f8f311SMark Adams #include <amgx_c.h>
13e6f8f311SMark Adams #include <limits>
14e6f8f311SMark Adams #include <vector>
15e6f8f311SMark Adams #include <algorithm>
16e6f8f311SMark Adams #include <map>
17e6f8f311SMark Adams #include <numeric>
18e6f8f311SMark Adams
199371c9d4SSatish Balay enum class AmgXSmoother {
209371c9d4SSatish Balay PCG,
219371c9d4SSatish Balay PCGF,
229371c9d4SSatish Balay PBiCGStab,
239371c9d4SSatish Balay GMRES,
249371c9d4SSatish Balay FGMRES,
259371c9d4SSatish Balay JacobiL1,
269371c9d4SSatish Balay BlockJacobi,
279371c9d4SSatish Balay GS,
289371c9d4SSatish Balay MulticolorGS,
299371c9d4SSatish Balay MulticolorILU,
309371c9d4SSatish Balay MulticolorDILU,
319371c9d4SSatish Balay ChebyshevPoly,
329371c9d4SSatish Balay NoSolver
339371c9d4SSatish Balay };
349371c9d4SSatish Balay enum class AmgXAMGMethod {
359371c9d4SSatish Balay Classical,
369371c9d4SSatish Balay Aggregation
379371c9d4SSatish Balay };
389371c9d4SSatish Balay enum class AmgXSelector {
399371c9d4SSatish Balay Size2,
409371c9d4SSatish Balay Size4,
419371c9d4SSatish Balay Size8,
429371c9d4SSatish Balay MultiPairwise,
439371c9d4SSatish Balay PMIS,
449371c9d4SSatish Balay HMIS
459371c9d4SSatish Balay };
469371c9d4SSatish Balay enum class AmgXCoarseSolver {
479371c9d4SSatish Balay DenseLU,
489371c9d4SSatish Balay NoSolver
499371c9d4SSatish Balay };
509371c9d4SSatish Balay enum class AmgXAMGCycle {
519371c9d4SSatish Balay V,
529371c9d4SSatish Balay W,
539371c9d4SSatish Balay F,
549371c9d4SSatish Balay CG,
559371c9d4SSatish Balay CGF
569371c9d4SSatish Balay };
57e6f8f311SMark Adams
589371c9d4SSatish Balay struct AmgXControlMap {
59e6f8f311SMark Adams static const std::map<std::string, AmgXAMGMethod> AMGMethods;
60e6f8f311SMark Adams static const std::map<std::string, AmgXSmoother> Smoothers;
61e6f8f311SMark Adams static const std::map<std::string, AmgXSelector> Selectors;
62e6f8f311SMark Adams static const std::map<std::string, AmgXCoarseSolver> CoarseSolvers;
63e6f8f311SMark Adams static const std::map<std::string, AmgXAMGCycle> AMGCycles;
64e6f8f311SMark Adams };
65e6f8f311SMark Adams
66e6f8f311SMark Adams const std::map<std::string, AmgXAMGMethod> AmgXControlMap::AMGMethods = {
67e6f8f311SMark Adams {"CLASSICAL", AmgXAMGMethod::Classical },
68e6f8f311SMark Adams {"AGGREGATION", AmgXAMGMethod::Aggregation}
69e6f8f311SMark Adams };
70e6f8f311SMark Adams
71e6f8f311SMark Adams const std::map<std::string, AmgXSmoother> AmgXControlMap::Smoothers = {
72e6f8f311SMark Adams {"PCG", AmgXSmoother::PCG },
73e6f8f311SMark Adams {"PCGF", AmgXSmoother::PCGF },
74e6f8f311SMark Adams {"PBICGSTAB", AmgXSmoother::PBiCGStab },
75e6f8f311SMark Adams {"GMRES", AmgXSmoother::GMRES },
76e6f8f311SMark Adams {"FGMRES", AmgXSmoother::FGMRES },
77e6f8f311SMark Adams {"JACOBI_L1", AmgXSmoother::JacobiL1 },
78e6f8f311SMark Adams {"BLOCK_JACOBI", AmgXSmoother::BlockJacobi },
79e6f8f311SMark Adams {"GS", AmgXSmoother::GS },
80e6f8f311SMark Adams {"MULTICOLOR_GS", AmgXSmoother::MulticolorGS },
81e6f8f311SMark Adams {"MULTICOLOR_ILU", AmgXSmoother::MulticolorILU },
82e6f8f311SMark Adams {"MULTICOLOR_DILU", AmgXSmoother::MulticolorDILU},
83e6f8f311SMark Adams {"CHEBYSHEV_POLY", AmgXSmoother::ChebyshevPoly },
84e6f8f311SMark Adams {"NOSOLVER", AmgXSmoother::NoSolver }
85e6f8f311SMark Adams };
86e6f8f311SMark Adams
87e6f8f311SMark Adams const std::map<std::string, AmgXSelector> AmgXControlMap::Selectors = {
88e6f8f311SMark Adams {"SIZE_2", AmgXSelector::Size2 },
89e6f8f311SMark Adams {"SIZE_4", AmgXSelector::Size4 },
90e6f8f311SMark Adams {"SIZE_8", AmgXSelector::Size8 },
91e6f8f311SMark Adams {"MULTI_PAIRWISE", AmgXSelector::MultiPairwise},
92e6f8f311SMark Adams {"PMIS", AmgXSelector::PMIS },
93e6f8f311SMark Adams {"HMIS", AmgXSelector::HMIS }
94e6f8f311SMark Adams };
95e6f8f311SMark Adams
96e6f8f311SMark Adams const std::map<std::string, AmgXCoarseSolver> AmgXControlMap::CoarseSolvers = {
97e6f8f311SMark Adams {"DENSE_LU_SOLVER", AmgXCoarseSolver::DenseLU },
98e6f8f311SMark Adams {"NOSOLVER", AmgXCoarseSolver::NoSolver}
99e6f8f311SMark Adams };
100e6f8f311SMark Adams
101e6f8f311SMark Adams const std::map<std::string, AmgXAMGCycle> AmgXControlMap::AMGCycles = {
102e6f8f311SMark Adams {"V", AmgXAMGCycle::V },
103e6f8f311SMark Adams {"W", AmgXAMGCycle::W },
104e6f8f311SMark Adams {"F", AmgXAMGCycle::F },
105e6f8f311SMark Adams {"CG", AmgXAMGCycle::CG },
106e6f8f311SMark Adams {"CGF", AmgXAMGCycle::CGF}
107e6f8f311SMark Adams };
108e6f8f311SMark Adams
109e6f8f311SMark Adams /*
110e6f8f311SMark Adams Private context (data structure) for the AMGX preconditioner.
111e6f8f311SMark Adams */
112e6f8f311SMark Adams struct PC_AMGX {
113e6f8f311SMark Adams AMGX_solver_handle solver;
114e6f8f311SMark Adams AMGX_config_handle cfg;
115e6f8f311SMark Adams AMGX_resources_handle rsrc;
116e6f8f311SMark Adams bool solve_state_init;
117e6f8f311SMark Adams bool rsrc_init;
118e6f8f311SMark Adams PetscBool verbose;
119e6f8f311SMark Adams
120e6f8f311SMark Adams AMGX_matrix_handle A;
121e6f8f311SMark Adams AMGX_vector_handle sol;
122e6f8f311SMark Adams AMGX_vector_handle rhs;
123e6f8f311SMark Adams
124e6f8f311SMark Adams MPI_Comm comm;
125e6f8f311SMark Adams PetscMPIInt rank = 0;
126e6f8f311SMark Adams PetscMPIInt nranks = 0;
127e6f8f311SMark Adams int devID = 0;
128e6f8f311SMark Adams
129e6f8f311SMark Adams void *lib_handle = 0;
130e6f8f311SMark Adams std::string cfg_contents;
131e6f8f311SMark Adams
132e6f8f311SMark Adams // Cached state for re-setup
133e6f8f311SMark Adams PetscInt nnz;
134e6f8f311SMark Adams PetscInt nLocalRows;
135e6f8f311SMark Adams PetscInt nGlobalRows;
136e6f8f311SMark Adams PetscInt bSize;
137e6f8f311SMark Adams Mat localA;
138e6f8f311SMark Adams const PetscScalar *values;
139e6f8f311SMark Adams
140e6f8f311SMark Adams // AMG Control parameters
141e6f8f311SMark Adams AmgXSmoother smoother;
142e6f8f311SMark Adams AmgXAMGMethod amg_method;
143e6f8f311SMark Adams AmgXSelector selector;
144e6f8f311SMark Adams AmgXCoarseSolver coarse_solver;
145e6f8f311SMark Adams AmgXAMGCycle amg_cycle;
146e6f8f311SMark Adams PetscInt presweeps;
147e6f8f311SMark Adams PetscInt postsweeps;
148e6f8f311SMark Adams PetscInt max_levels;
149e6f8f311SMark Adams PetscInt aggressive_levels;
150e6f8f311SMark Adams PetscInt dense_lu_num_rows;
151e6f8f311SMark Adams PetscScalar strength_threshold;
152e6f8f311SMark Adams PetscBool print_grid_stats;
153e6f8f311SMark Adams PetscBool exact_coarse_solve;
154e6f8f311SMark Adams
155e6f8f311SMark Adams // Smoother control parameters
156e6f8f311SMark Adams PetscScalar jacobi_relaxation_factor;
157e6f8f311SMark Adams PetscScalar gs_symmetric;
158e6f8f311SMark Adams };
159e6f8f311SMark Adams
160e6f8f311SMark Adams static PetscInt s_count = 0;
161e6f8f311SMark Adams
162e6f8f311SMark Adams // Buffer of messages from AmgX
163e6f8f311SMark Adams // Currently necessary hack before we adapt AmgX to print from single rank only
164e6f8f311SMark Adams static std::string amgx_output{};
165e6f8f311SMark Adams
166e6f8f311SMark Adams // A print callback that allows AmgX to return status messages
print_callback(const char * msg,int length)167d71ae5a4SJacob Faibussowitsch static void print_callback(const char *msg, int length)
168d71ae5a4SJacob Faibussowitsch {
169e6f8f311SMark Adams amgx_output.append(msg);
170e6f8f311SMark Adams }
171e6f8f311SMark Adams
172e6f8f311SMark Adams // Outputs messages from the AmgX message buffer and clears it
amgx_output_messages(PC_AMGX * amgx)17366976f2fSJacob Faibussowitsch static PetscErrorCode amgx_output_messages(PC_AMGX *amgx)
174d71ae5a4SJacob Faibussowitsch {
175e6f8f311SMark Adams PetscFunctionBegin;
176e6f8f311SMark Adams // If AmgX output is enabled and we have a message, output it
177e6f8f311SMark Adams if (amgx->verbose && !amgx_output.empty()) {
178e6f8f311SMark Adams // Only a single rank to output the AmgX messages
179e6f8f311SMark Adams PetscCall(PetscPrintf(amgx->comm, "AMGX: %s", amgx_output.c_str()));
180e6f8f311SMark Adams
181e6f8f311SMark Adams // Note that all ranks clear their received output
182e6f8f311SMark Adams amgx_output.clear();
183e6f8f311SMark Adams }
1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
185e6f8f311SMark Adams }
186e6f8f311SMark Adams
187e6f8f311SMark Adams // XXX Need to add call in AmgX API that gracefully destroys everything
188e6f8f311SMark Adams // without abort etc.
1899371c9d4SSatish Balay #define PetscCallAmgX(rc) \
1909371c9d4SSatish Balay do { \
191e6f8f311SMark Adams AMGX_RC err = (rc); \
192e6f8f311SMark Adams char msg[4096]; \
193e6f8f311SMark Adams switch (err) { \
194d71ae5a4SJacob Faibussowitsch case AMGX_RC_OK: \
195d71ae5a4SJacob Faibussowitsch break; \
196d71ae5a4SJacob Faibussowitsch default: \
197d71ae5a4SJacob Faibussowitsch AMGX_get_error_string(err, msg, 4096); \
198d71ae5a4SJacob Faibussowitsch SETERRQ(amgx->comm, PETSC_ERR_LIB, "%s", msg); \
199e6f8f311SMark Adams } \
200e6f8f311SMark Adams } while (0)
201e6f8f311SMark Adams
202e6f8f311SMark Adams /*
203e6f8f311SMark Adams PCSetUp_AMGX - Prepares for the use of the AmgX preconditioner
204e6f8f311SMark Adams by setting data structures and options.
205e6f8f311SMark Adams
206e6f8f311SMark Adams Input Parameter:
207e6f8f311SMark Adams . pc - the preconditioner context
208e6f8f311SMark Adams
209e6f8f311SMark Adams Application Interface Routine: PCSetUp()
210e6f8f311SMark Adams
211f1580f4eSBarry Smith Note:
212e6f8f311SMark Adams The interface routine PCSetUp() is not usually called directly by
213e6f8f311SMark Adams the user, but instead is called by PCApply() if necessary.
214e6f8f311SMark Adams */
PCSetUp_AMGX(PC pc)215d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_AMGX(PC pc)
216d71ae5a4SJacob Faibussowitsch {
217e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
218e6f8f311SMark Adams Mat Pmat = pc->pmat;
219e6f8f311SMark Adams PetscBool is_dev_ptrs;
220e6f8f311SMark Adams
221e6f8f311SMark Adams PetscFunctionBegin;
222e6f8f311SMark Adams PetscCall(PetscObjectTypeCompareAny((PetscObject)Pmat, &is_dev_ptrs, MATAIJCUSPARSE, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, ""));
223e6f8f311SMark Adams
224e6f8f311SMark Adams // At the present time, an AmgX matrix is a sequential matrix
225e6f8f311SMark Adams // Non-sequential/MPI matrices must be adapted to extract the local matrix
226e6f8f311SMark Adams bool partial_setup_allowed = (pc->setupcalled && pc->flag != DIFFERENT_NONZERO_PATTERN);
227e6f8f311SMark Adams if (amgx->nranks > 1) {
228e6f8f311SMark Adams if (partial_setup_allowed) {
229e6f8f311SMark Adams PetscCall(MatMPIAIJGetLocalMat(Pmat, MAT_REUSE_MATRIX, &amgx->localA));
230e6f8f311SMark Adams } else {
231e6f8f311SMark Adams PetscCall(MatMPIAIJGetLocalMat(Pmat, MAT_INITIAL_MATRIX, &amgx->localA));
232e6f8f311SMark Adams }
233e6f8f311SMark Adams
23448a46eb9SPierre Jolivet if (is_dev_ptrs) PetscCall(MatConvert(amgx->localA, MATSEQAIJCUSPARSE, MAT_INPLACE_MATRIX, &amgx->localA));
235e6f8f311SMark Adams } else {
236e6f8f311SMark Adams amgx->localA = Pmat;
237e6f8f311SMark Adams }
238e6f8f311SMark Adams
239e6f8f311SMark Adams if (is_dev_ptrs) {
240e6f8f311SMark Adams PetscCall(MatSeqAIJCUSPARSEGetArrayRead(amgx->localA, &amgx->values));
241e6f8f311SMark Adams } else {
242e6f8f311SMark Adams PetscCall(MatSeqAIJGetArrayRead(amgx->localA, &amgx->values));
243e6f8f311SMark Adams }
244e6f8f311SMark Adams
245e6f8f311SMark Adams if (!partial_setup_allowed) {
246e6f8f311SMark Adams // Initialise resources and matrices
247e6f8f311SMark Adams if (!amgx->rsrc_init) {
248e6f8f311SMark Adams // Read configuration file
249e6f8f311SMark Adams PetscCallAmgX(AMGX_config_create(&amgx->cfg, amgx->cfg_contents.c_str()));
250e6f8f311SMark Adams PetscCallAmgX(AMGX_resources_create(&amgx->rsrc, amgx->cfg, &amgx->comm, 1, &amgx->devID));
251e6f8f311SMark Adams amgx->rsrc_init = true;
252e6f8f311SMark Adams }
253e6f8f311SMark Adams
254e6f8f311SMark Adams PetscCheck(!amgx->solve_state_init, amgx->comm, PETSC_ERR_PLIB, "AmgX solve state initialisation already called.");
255e6f8f311SMark Adams PetscCallAmgX(AMGX_matrix_create(&amgx->A, amgx->rsrc, AMGX_mode_dDDI));
256e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_create(&amgx->sol, amgx->rsrc, AMGX_mode_dDDI));
257e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_create(&amgx->rhs, amgx->rsrc, AMGX_mode_dDDI));
258e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_create(&amgx->solver, amgx->rsrc, AMGX_mode_dDDI, amgx->cfg));
259e6f8f311SMark Adams amgx->solve_state_init = true;
260e6f8f311SMark Adams
261e6f8f311SMark Adams // Extract the CSR data
262e6f8f311SMark Adams PetscBool done;
263e6f8f311SMark Adams const PetscInt *colIndices;
264e6f8f311SMark Adams const PetscInt *rowOffsets;
265e6f8f311SMark Adams PetscCall(MatGetRowIJ(amgx->localA, 0, PETSC_FALSE, PETSC_FALSE, &amgx->nLocalRows, &rowOffsets, &colIndices, &done));
266e6f8f311SMark Adams PetscCheck(done, amgx->comm, PETSC_ERR_PLIB, "MatGetRowIJ was not successful");
267e6f8f311SMark Adams PetscCheck(amgx->nLocalRows < std::numeric_limits<int>::max(), PETSC_COMM_SELF, PETSC_ERR_PLIB, "AmgX restricted to int local rows but nLocalRows = %" PetscInt_FMT " > max<int>", amgx->nLocalRows);
268e6f8f311SMark Adams
269e6f8f311SMark Adams if (is_dev_ptrs) {
270e6f8f311SMark Adams PetscCallCUDA(cudaMemcpy(&amgx->nnz, &rowOffsets[amgx->nLocalRows], sizeof(int), cudaMemcpyDefault));
271e6f8f311SMark Adams } else {
272e6f8f311SMark Adams amgx->nnz = rowOffsets[amgx->nLocalRows];
273e6f8f311SMark Adams }
274e6f8f311SMark Adams
275e6f8f311SMark Adams PetscCheck(amgx->nnz < std::numeric_limits<int>::max(), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support for 64-bit integer nnz not yet implemented, nnz = %" PetscInt_FMT ".", amgx->nnz);
276e6f8f311SMark Adams
277e6f8f311SMark Adams // Allocate space for some partition offsets
278e6f8f311SMark Adams std::vector<PetscInt> partitionOffsets(amgx->nranks + 1);
279e6f8f311SMark Adams
280e6f8f311SMark Adams // Fetch the number of local rows per rank
281e6f8f311SMark Adams partitionOffsets[0] = 0; /* could use PetscLayoutGetRanges */
282e6f8f311SMark Adams PetscCallMPI(MPI_Allgather(&amgx->nLocalRows, 1, MPIU_INT, partitionOffsets.data() + 1, 1, MPIU_INT, amgx->comm));
283e6f8f311SMark Adams std::partial_sum(partitionOffsets.begin(), partitionOffsets.end(), partitionOffsets.begin());
284e6f8f311SMark Adams
285e6f8f311SMark Adams // Fetch the number of global rows
286e6f8f311SMark Adams amgx->nGlobalRows = partitionOffsets[amgx->nranks];
287e6f8f311SMark Adams
288e6f8f311SMark Adams PetscCall(MatGetBlockSize(Pmat, &amgx->bSize));
289e6f8f311SMark Adams
290e6f8f311SMark Adams // XXX Currently constrained to 32-bit indices, to be changed in the future
291e6f8f311SMark Adams // Create the distribution and upload the matrix data
292e6f8f311SMark Adams AMGX_distribution_handle dist;
293e6f8f311SMark Adams PetscCallAmgX(AMGX_distribution_create(&dist, amgx->cfg));
294e6f8f311SMark Adams PetscCallAmgX(AMGX_distribution_set_32bit_colindices(dist, true));
295e6f8f311SMark Adams PetscCallAmgX(AMGX_distribution_set_partition_data(dist, AMGX_DIST_PARTITION_OFFSETS, partitionOffsets.data()));
296e6f8f311SMark Adams PetscCallAmgX(AMGX_matrix_upload_distributed(amgx->A, amgx->nGlobalRows, (int)amgx->nLocalRows, (int)amgx->nnz, amgx->bSize, amgx->bSize, rowOffsets, colIndices, amgx->values, NULL, dist));
297e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_setup(amgx->solver, amgx->A));
298e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_bind(amgx->sol, amgx->A));
299e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_bind(amgx->rhs, amgx->A));
300e6f8f311SMark Adams
301e6f8f311SMark Adams PetscInt nlr = 0;
302e6f8f311SMark Adams PetscCall(MatRestoreRowIJ(amgx->localA, 0, PETSC_FALSE, PETSC_FALSE, &nlr, &rowOffsets, &colIndices, &done));
303e6f8f311SMark Adams } else {
304e6f8f311SMark Adams // The fast path for if the sparsity pattern persists
305e6f8f311SMark Adams PetscCallAmgX(AMGX_matrix_replace_coefficients(amgx->A, amgx->nLocalRows, amgx->nnz, amgx->values, NULL));
306e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_resetup(amgx->solver, amgx->A));
307e6f8f311SMark Adams }
308e6f8f311SMark Adams
309e6f8f311SMark Adams if (is_dev_ptrs) {
310e6f8f311SMark Adams PetscCall(MatSeqAIJCUSPARSERestoreArrayRead(amgx->localA, &amgx->values));
311e6f8f311SMark Adams } else {
312e6f8f311SMark Adams PetscCall(MatSeqAIJRestoreArrayRead(amgx->localA, &amgx->values));
313e6f8f311SMark Adams }
3143ba16761SJacob Faibussowitsch PetscCall(amgx_output_messages(amgx));
3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
316e6f8f311SMark Adams }
317e6f8f311SMark Adams
318e6f8f311SMark Adams /*
319e6f8f311SMark Adams PCApply_AMGX - Applies the AmgX preconditioner to a vector.
320e6f8f311SMark Adams
321e6f8f311SMark Adams Input Parameters:
322e6f8f311SMark Adams . pc - the preconditioner context
323e6f8f311SMark Adams . b - rhs vector
324e6f8f311SMark Adams
325e6f8f311SMark Adams Output Parameter:
326e6f8f311SMark Adams . x - solution vector
327e6f8f311SMark Adams
328e6f8f311SMark Adams Application Interface Routine: PCApply()
329e6f8f311SMark Adams */
PCApply_AMGX(PC pc,Vec b,Vec x)330d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_AMGX(PC pc, Vec b, Vec x)
331d71ae5a4SJacob Faibussowitsch {
332e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
333e6f8f311SMark Adams PetscScalar *x_;
334e6f8f311SMark Adams const PetscScalar *b_;
335e6f8f311SMark Adams PetscBool is_dev_ptrs;
336e6f8f311SMark Adams
337e6f8f311SMark Adams PetscFunctionBegin;
338e6f8f311SMark Adams PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &is_dev_ptrs, VECCUDA, VECMPICUDA, VECSEQCUDA, ""));
339e6f8f311SMark Adams
340e6f8f311SMark Adams if (is_dev_ptrs) {
341e6f8f311SMark Adams PetscCall(VecCUDAGetArrayWrite(x, &x_));
342e6f8f311SMark Adams PetscCall(VecCUDAGetArrayRead(b, &b_));
343e6f8f311SMark Adams } else {
344e6f8f311SMark Adams PetscCall(VecGetArrayWrite(x, &x_));
345e6f8f311SMark Adams PetscCall(VecGetArrayRead(b, &b_));
346e6f8f311SMark Adams }
347e6f8f311SMark Adams
348e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_upload(amgx->sol, amgx->nLocalRows, 1, x_));
349e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_upload(amgx->rhs, amgx->nLocalRows, 1, b_));
350e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_solve_with_0_initial_guess(amgx->solver, amgx->rhs, amgx->sol));
351e6f8f311SMark Adams
352e6f8f311SMark Adams AMGX_SOLVE_STATUS status;
353e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_get_status(amgx->solver, &status));
354e6f8f311SMark Adams PetscCall(PCSetErrorIfFailure(pc, static_cast<PetscBool>(status == AMGX_SOLVE_FAILED)));
355e6f8f311SMark Adams PetscCheck(status != AMGX_SOLVE_FAILED, amgx->comm, PETSC_ERR_CONV_FAILED, "AmgX solver failed to solve the system! The error code is %d.", status);
356e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_download(amgx->sol, x_));
357e6f8f311SMark Adams
358e6f8f311SMark Adams if (is_dev_ptrs) {
359e6f8f311SMark Adams PetscCall(VecCUDARestoreArrayWrite(x, &x_));
360e6f8f311SMark Adams PetscCall(VecCUDARestoreArrayRead(b, &b_));
361e6f8f311SMark Adams } else {
362e6f8f311SMark Adams PetscCall(VecRestoreArrayWrite(x, &x_));
363e6f8f311SMark Adams PetscCall(VecRestoreArrayRead(b, &b_));
364e6f8f311SMark Adams }
3653ba16761SJacob Faibussowitsch PetscCall(amgx_output_messages(amgx));
3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
367e6f8f311SMark Adams }
368e6f8f311SMark Adams
PCReset_AMGX(PC pc)369d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_AMGX(PC pc)
370d71ae5a4SJacob Faibussowitsch {
371e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
372e6f8f311SMark Adams
373e6f8f311SMark Adams PetscFunctionBegin;
374e6f8f311SMark Adams if (amgx->solve_state_init) {
375e6f8f311SMark Adams PetscCallAmgX(AMGX_solver_destroy(amgx->solver));
376e6f8f311SMark Adams PetscCallAmgX(AMGX_matrix_destroy(amgx->A));
377e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_destroy(amgx->sol));
378e6f8f311SMark Adams PetscCallAmgX(AMGX_vector_destroy(amgx->rhs));
379e6f8f311SMark Adams if (amgx->nranks > 1) PetscCall(MatDestroy(&amgx->localA));
3803ba16761SJacob Faibussowitsch PetscCall(amgx_output_messages(amgx));
381e6f8f311SMark Adams amgx->solve_state_init = false;
382e6f8f311SMark Adams }
3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
384e6f8f311SMark Adams }
385e6f8f311SMark Adams
386e6f8f311SMark Adams /*
387e6f8f311SMark Adams PCDestroy_AMGX - Destroys the private context for the AmgX preconditioner
388e6f8f311SMark Adams that was created with PCCreate_AMGX().
389e6f8f311SMark Adams
390e6f8f311SMark Adams Input Parameter:
391e6f8f311SMark Adams . pc - the preconditioner context
392e6f8f311SMark Adams
393e6f8f311SMark Adams Application Interface Routine: PCDestroy()
394e6f8f311SMark Adams */
PCDestroy_AMGX(PC pc)395d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_AMGX(PC pc)
396d71ae5a4SJacob Faibussowitsch {
397e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
398e6f8f311SMark Adams
399e6f8f311SMark Adams PetscFunctionBegin;
400e6f8f311SMark Adams /* decrease the number of instances, only the last instance need to destroy resource and finalizing AmgX */
401e6f8f311SMark Adams if (s_count == 1) {
402e6f8f311SMark Adams /* can put this in a PCAMGXInitializePackage method */
403e6f8f311SMark Adams PetscCheck(amgx->rsrc != nullptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "s_rsrc == NULL");
404e6f8f311SMark Adams PetscCallAmgX(AMGX_resources_destroy(amgx->rsrc));
405e6f8f311SMark Adams /* destroy config (need to use AMGX_SAFE_CALL after this point) */
406e6f8f311SMark Adams PetscCallAmgX(AMGX_config_destroy(amgx->cfg));
407e6f8f311SMark Adams PetscCallAmgX(AMGX_finalize_plugins());
408e6f8f311SMark Adams PetscCallAmgX(AMGX_finalize());
409e6f8f311SMark Adams PetscCallMPI(MPI_Comm_free(&amgx->comm));
410e6f8f311SMark Adams } else {
411e6f8f311SMark Adams PetscCallAmgX(AMGX_config_destroy(amgx->cfg));
412e6f8f311SMark Adams }
413e6f8f311SMark Adams s_count -= 1;
414e6f8f311SMark Adams PetscCall(PetscFree(amgx));
4153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
416e6f8f311SMark Adams }
417e6f8f311SMark Adams
418e6f8f311SMark Adams template <class T>
map_reverse_lookup(const std::map<std::string,T> & map,const T & key)419d71ae5a4SJacob Faibussowitsch std::string map_reverse_lookup(const std::map<std::string, T> &map, const T &key)
420d71ae5a4SJacob Faibussowitsch {
421e6f8f311SMark Adams for (auto const &m : map) {
422ad540459SPierre Jolivet if (m.second == key) return m.first;
423e6f8f311SMark Adams }
424e6f8f311SMark Adams return "";
425e6f8f311SMark Adams }
426e6f8f311SMark Adams
PCSetFromOptions_AMGX(PC pc,PetscOptionItems PetscOptionsObject)427ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_AMGX(PC pc, PetscOptionItems PetscOptionsObject)
428d71ae5a4SJacob Faibussowitsch {
429e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
430e6f8f311SMark Adams constexpr int MAX_PARAM_LEN = 128;
431e6f8f311SMark Adams char option[MAX_PARAM_LEN];
432e6f8f311SMark Adams
433e6f8f311SMark Adams PetscFunctionBegin;
434e6f8f311SMark Adams PetscOptionsHeadBegin(PetscOptionsObject, "AmgX options");
435e6f8f311SMark Adams amgx->cfg_contents = "config_version=2,";
436e6f8f311SMark Adams amgx->cfg_contents += "determinism_flag=1,";
437e6f8f311SMark Adams
438e6f8f311SMark Adams // Set exact coarse solve
439e6f8f311SMark Adams PetscCall(PetscOptionsBool("-pc_amgx_exact_coarse_solve", "AmgX AMG Exact Coarse Solve", "", amgx->exact_coarse_solve, &amgx->exact_coarse_solve, NULL));
440ad540459SPierre Jolivet if (amgx->exact_coarse_solve) amgx->cfg_contents += "exact_coarse_solve=1,";
441e6f8f311SMark Adams
442e6f8f311SMark Adams amgx->cfg_contents += "solver(amg)=AMG,";
443e6f8f311SMark Adams
444e6f8f311SMark Adams // Set method
445e6f8f311SMark Adams std::string def_amg_method = map_reverse_lookup(AmgXControlMap::AMGMethods, amgx->amg_method);
446c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, def_amg_method.c_str(), sizeof(option)));
447e6f8f311SMark Adams PetscCall(PetscOptionsString("-pc_amgx_amg_method", "AmgX AMG Method", "", option, option, MAX_PARAM_LEN, NULL));
448e6f8f311SMark Adams PetscCheck(AmgXControlMap::AMGMethods.count(option) == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "AMG Method %s not registered for AmgX.", option);
449e6f8f311SMark Adams amgx->amg_method = AmgXControlMap::AMGMethods.at(option);
450e6f8f311SMark Adams amgx->cfg_contents += "amg:algorithm=" + std::string(option) + ",";
451e6f8f311SMark Adams
452e6f8f311SMark Adams // Set cycle
453e6f8f311SMark Adams std::string def_amg_cycle = map_reverse_lookup(AmgXControlMap::AMGCycles, amgx->amg_cycle);
454c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, def_amg_cycle.c_str(), sizeof(option)));
455e6f8f311SMark Adams PetscCall(PetscOptionsString("-pc_amgx_amg_cycle", "AmgX AMG Cycle", "", option, option, MAX_PARAM_LEN, NULL));
456e6f8f311SMark Adams PetscCheck(AmgXControlMap::AMGCycles.count(option) == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "AMG Cycle %s not registered for AmgX.", option);
457e6f8f311SMark Adams amgx->amg_cycle = AmgXControlMap::AMGCycles.at(option);
458e6f8f311SMark Adams amgx->cfg_contents += "amg:cycle=" + std::string(option) + ",";
459e6f8f311SMark Adams
460e6f8f311SMark Adams // Set smoother
461e6f8f311SMark Adams std::string def_smoother = map_reverse_lookup(AmgXControlMap::Smoothers, amgx->smoother);
462c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, def_smoother.c_str(), sizeof(option)));
463e6f8f311SMark Adams PetscCall(PetscOptionsString("-pc_amgx_smoother", "AmgX Smoother", "", option, option, MAX_PARAM_LEN, NULL));
464e6f8f311SMark Adams PetscCheck(AmgXControlMap::Smoothers.count(option) == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Smoother %s not registered for AmgX.", option);
465e6f8f311SMark Adams amgx->smoother = AmgXControlMap::Smoothers.at(option);
466e6f8f311SMark Adams amgx->cfg_contents += "amg:smoother(smooth)=" + std::string(option) + ",";
467e6f8f311SMark Adams
468e6f8f311SMark Adams if (amgx->smoother == AmgXSmoother::JacobiL1 || amgx->smoother == AmgXSmoother::BlockJacobi) {
469e6f8f311SMark Adams PetscCall(PetscOptionsScalar("-pc_amgx_jacobi_relaxation_factor", "AmgX AMG Jacobi Relaxation Factor", "", amgx->jacobi_relaxation_factor, &amgx->jacobi_relaxation_factor, NULL));
470e6f8f311SMark Adams amgx->cfg_contents += "smooth:relaxation_factor=" + std::to_string(amgx->jacobi_relaxation_factor) + ",";
471e6f8f311SMark Adams } else if (amgx->smoother == AmgXSmoother::GS || amgx->smoother == AmgXSmoother::MulticolorGS) {
472e6f8f311SMark Adams PetscCall(PetscOptionsScalar("-pc_amgx_gs_symmetric", "AmgX AMG Gauss Seidel Symmetric", "", amgx->gs_symmetric, &amgx->gs_symmetric, NULL));
473e6f8f311SMark Adams amgx->cfg_contents += "smooth:symmetric_GS=" + std::to_string(amgx->gs_symmetric) + ",";
474e6f8f311SMark Adams }
475e6f8f311SMark Adams
476e6f8f311SMark Adams // Set selector
477e6f8f311SMark Adams std::string def_selector = map_reverse_lookup(AmgXControlMap::Selectors, amgx->selector);
478c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, def_selector.c_str(), sizeof(option)));
479e6f8f311SMark Adams PetscCall(PetscOptionsString("-pc_amgx_selector", "AmgX Selector", "", option, option, MAX_PARAM_LEN, NULL));
480e6f8f311SMark Adams PetscCheck(AmgXControlMap::Selectors.count(option) == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Selector %s not registered for AmgX.", option);
481e6f8f311SMark Adams
482e6f8f311SMark Adams // Double check that the user has selected an appropriate selector for the AMG method
483e6f8f311SMark Adams if (amgx->amg_method == AmgXAMGMethod::Classical) {
484e6f8f311SMark Adams PetscCheck(amgx->selector == AmgXSelector::PMIS || amgx->selector == AmgXSelector::HMIS, amgx->comm, PETSC_ERR_PLIB, "Chosen selector is not used for AmgX Classical AMG: selector=%s", option);
485e6f8f311SMark Adams amgx->cfg_contents += "amg:interpolator=D2,";
486e6f8f311SMark Adams } else if (amgx->amg_method == AmgXAMGMethod::Aggregation) {
487e6f8f311SMark Adams PetscCheck(amgx->selector == AmgXSelector::Size2 || amgx->selector == AmgXSelector::Size4 || amgx->selector == AmgXSelector::Size8 || amgx->selector == AmgXSelector::MultiPairwise, amgx->comm, PETSC_ERR_PLIB, "Chosen selector is not used for AmgX Aggregation AMG");
488e6f8f311SMark Adams }
489e6f8f311SMark Adams amgx->selector = AmgXControlMap::Selectors.at(option);
490e6f8f311SMark Adams amgx->cfg_contents += "amg:selector=" + std::string(option) + ",";
491e6f8f311SMark Adams
492e6f8f311SMark Adams // Set presweeps
493e6f8f311SMark Adams PetscCall(PetscOptionsInt("-pc_amgx_presweeps", "AmgX AMG Presweep Count", "", amgx->presweeps, &amgx->presweeps, NULL));
494e6f8f311SMark Adams amgx->cfg_contents += "amg:presweeps=" + std::to_string(amgx->presweeps) + ",";
495e6f8f311SMark Adams
496e6f8f311SMark Adams // Set postsweeps
497e6f8f311SMark Adams PetscCall(PetscOptionsInt("-pc_amgx_postsweeps", "AmgX AMG Postsweep Count", "", amgx->postsweeps, &amgx->postsweeps, NULL));
498e6f8f311SMark Adams amgx->cfg_contents += "amg:postsweeps=" + std::to_string(amgx->postsweeps) + ",";
499e6f8f311SMark Adams
500e6f8f311SMark Adams // Set max levels
501e6f8f311SMark Adams PetscCall(PetscOptionsInt("-pc_amgx_max_levels", "AmgX AMG Max Level Count", "", amgx->max_levels, &amgx->max_levels, NULL));
5020bdf8ccaSPierre Jolivet amgx->cfg_contents += "amg:max_levels=" + std::to_string(amgx->max_levels) + ",";
503e6f8f311SMark Adams
504e6f8f311SMark Adams // Set dense LU num rows
505e6f8f311SMark Adams PetscCall(PetscOptionsInt("-pc_amgx_dense_lu_num_rows", "AmgX Dense LU Number of Rows", "", amgx->dense_lu_num_rows, &amgx->dense_lu_num_rows, NULL));
506e6f8f311SMark Adams amgx->cfg_contents += "amg:dense_lu_num_rows=" + std::to_string(amgx->dense_lu_num_rows) + ",";
507e6f8f311SMark Adams
508e6f8f311SMark Adams // Set strength threshold
509e6f8f311SMark Adams PetscCall(PetscOptionsScalar("-pc_amgx_strength_threshold", "AmgX AMG Strength Threshold", "", amgx->strength_threshold, &amgx->strength_threshold, NULL));
510e6f8f311SMark Adams amgx->cfg_contents += "amg:strength_threshold=" + std::to_string(amgx->strength_threshold) + ",";
511e6f8f311SMark Adams
512e6f8f311SMark Adams // Set aggressive_levels
513e6f8f311SMark Adams PetscCall(PetscOptionsInt("-pc_amgx_aggressive_levels", "AmgX AMG Presweep Count", "", amgx->aggressive_levels, &amgx->aggressive_levels, NULL));
514ad540459SPierre Jolivet if (amgx->aggressive_levels > 0) amgx->cfg_contents += "amg:aggressive_levels=" + std::to_string(amgx->aggressive_levels) + ",";
515e6f8f311SMark Adams
516e6f8f311SMark Adams // Set coarse solver
517e6f8f311SMark Adams std::string def_coarse_solver = map_reverse_lookup(AmgXControlMap::CoarseSolvers, amgx->coarse_solver);
518c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, def_coarse_solver.c_str(), sizeof(option)));
519e6f8f311SMark Adams PetscCall(PetscOptionsString("-pc_amgx_coarse_solver", "AmgX CoarseSolver", "", option, option, MAX_PARAM_LEN, NULL));
520e6f8f311SMark Adams PetscCheck(AmgXControlMap::CoarseSolvers.count(option) == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "CoarseSolver %s not registered for AmgX.", option);
521e6f8f311SMark Adams amgx->coarse_solver = AmgXControlMap::CoarseSolvers.at(option);
522e6f8f311SMark Adams amgx->cfg_contents += "amg:coarse_solver=" + std::string(option) + ",";
523e6f8f311SMark Adams
524e6f8f311SMark Adams // Set max iterations
525e6f8f311SMark Adams amgx->cfg_contents += "amg:max_iters=1,";
526e6f8f311SMark Adams
527e6f8f311SMark Adams // Set output control parameters
528e6f8f311SMark Adams PetscCall(PetscOptionsBool("-pc_amgx_print_grid_stats", "AmgX Print Grid Stats", "", amgx->print_grid_stats, &amgx->print_grid_stats, NULL));
529e6f8f311SMark Adams
530ad540459SPierre Jolivet if (amgx->print_grid_stats) amgx->cfg_contents += "amg:print_grid_stats=1,";
531e6f8f311SMark Adams amgx->cfg_contents += "amg:monitor_residual=0";
532e6f8f311SMark Adams
533e6f8f311SMark Adams // Set whether AmgX output will be seen
534e6f8f311SMark Adams PetscCall(PetscOptionsBool("-pc_amgx_verbose", "Enable output from AmgX", "", amgx->verbose, &amgx->verbose, NULL));
535e6f8f311SMark Adams PetscOptionsHeadEnd();
5363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
537e6f8f311SMark Adams }
538e6f8f311SMark Adams
PCView_AMGX(PC pc,PetscViewer viewer)539d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_AMGX(PC pc, PetscViewer viewer)
540d71ae5a4SJacob Faibussowitsch {
541e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
542*9f196a02SMartin Diehl PetscBool isascii;
543e6f8f311SMark Adams
544e6f8f311SMark Adams PetscFunctionBegin;
545*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
546*9f196a02SMartin Diehl if (isascii) {
547e6f8f311SMark Adams std::string output_cfg(amgx->cfg_contents);
548e6f8f311SMark Adams std::replace(output_cfg.begin(), output_cfg.end(), ',', '\n');
549e6f8f311SMark Adams PetscCall(PetscViewerASCIIPrintf(viewer, "\n%s\n", output_cfg.c_str()));
550e6f8f311SMark Adams }
5513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
552e6f8f311SMark Adams }
553e6f8f311SMark Adams
554e6f8f311SMark Adams /*MC
555e6f8f311SMark Adams PCAMGX - Interface to NVIDIA's AmgX algebraic multigrid
556e6f8f311SMark Adams
557e6f8f311SMark Adams Options Database Keys:
558e6f8f311SMark Adams + -pc_amgx_amg_method <CLASSICAL,AGGREGATION> - set the AMG algorithm to use
559e6f8f311SMark Adams . -pc_amgx_amg_cycle <V,W,F,CG> - set the AMG cycle type
560e6f8f311SMark Adams . -pc_amgx_smoother <PCG,PCGF,PBICGSTAB,GMRES,FGMRES,JACOBI_L1,BLOCK_JACOBI,GS,MULTICOLOR_GS,MULTICOLOR_ILU,MULTICOLOR_DILU,CHEBYSHEV_POLY,NOSOLVER> - set the AMG pre/post smoother
561e6f8f311SMark Adams . -pc_amgx_jacobi_relaxation_factor - set the relaxation factor for Jacobi smoothing
562e6f8f311SMark Adams . -pc_amgx_gs_symmetric - enforce symmetric Gauss-Seidel smoothing (only applies if GS smoothing is selected)
563e6f8f311SMark Adams . -pc_amgx_selector <SIZE_2,SIZE_4,SIZE_8,MULTI_PAIRWISE,PMIS,HMIS> - set the AMG coarse selector
564e6f8f311SMark Adams . -pc_amgx_presweeps - set the number of AMG pre-sweeps
565e6f8f311SMark Adams . -pc_amgx_postsweeps - set the number of AMG post-sweeps
566e6f8f311SMark Adams . -pc_amgx_max_levels - set the maximum number of levels in the AMG level hierarchy
567e6f8f311SMark Adams . -pc_amgx_strength_threshold - set the strength threshold for the AMG coarsening
568e6f8f311SMark Adams . -pc_amgx_aggressive_levels - set the number of levels (from the finest) that should apply aggressive coarsening
569e6f8f311SMark Adams . -pc_amgx_coarse_solver <DENSE_LU_SOLVER,NOSOLVER> - set the coarse solve
570e6f8f311SMark Adams . -pc_amgx_print_grid_stats - output the AMG grid hierarchy to stdout
571e6f8f311SMark Adams - -pc_amgx_verbose - enable AmgX output
572e6f8f311SMark Adams
573e6f8f311SMark Adams Level: intermediate
574e6f8f311SMark Adams
575f1580f4eSBarry Smith Note:
576f1580f4eSBarry Smith Implementation will accept host or device pointers, but good performance will require that the `KSP` is also GPU accelerated so that data is not frequently transferred between host and device.
577e6f8f311SMark Adams
578562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCHYPRE`, `PCMG`, `PCAmgXGetResources()`, `PCCreate()`, `PCSetType()`, `PCType` (for list of available types), `PC`
579e6f8f311SMark Adams M*/
580e6f8f311SMark Adams
PCCreate_AMGX(PC pc)581d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_AMGX(PC pc)
582d71ae5a4SJacob Faibussowitsch {
583e6f8f311SMark Adams PC_AMGX *amgx;
584e6f8f311SMark Adams
585e6f8f311SMark Adams PetscFunctionBegin;
5864dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&amgx));
587e6f8f311SMark Adams pc->ops->apply = PCApply_AMGX;
588e6f8f311SMark Adams pc->ops->setfromoptions = PCSetFromOptions_AMGX;
589e6f8f311SMark Adams pc->ops->setup = PCSetUp_AMGX;
590e6f8f311SMark Adams pc->ops->view = PCView_AMGX;
591e6f8f311SMark Adams pc->ops->destroy = PCDestroy_AMGX;
592e6f8f311SMark Adams pc->ops->reset = PCReset_AMGX;
593e6f8f311SMark Adams pc->data = (void *)amgx;
594e6f8f311SMark Adams
595e6f8f311SMark Adams // Set the defaults
596e6f8f311SMark Adams amgx->selector = AmgXSelector::PMIS;
597e6f8f311SMark Adams amgx->smoother = AmgXSmoother::BlockJacobi;
598e6f8f311SMark Adams amgx->amg_method = AmgXAMGMethod::Classical;
599e6f8f311SMark Adams amgx->coarse_solver = AmgXCoarseSolver::DenseLU;
600e6f8f311SMark Adams amgx->amg_cycle = AmgXAMGCycle::V;
601e6f8f311SMark Adams amgx->exact_coarse_solve = PETSC_TRUE;
602e6f8f311SMark Adams amgx->presweeps = 1;
603e6f8f311SMark Adams amgx->postsweeps = 1;
604e6f8f311SMark Adams amgx->max_levels = 100;
605e6f8f311SMark Adams amgx->strength_threshold = 0.5;
606e6f8f311SMark Adams amgx->aggressive_levels = 0;
607e6f8f311SMark Adams amgx->dense_lu_num_rows = 1;
608e6f8f311SMark Adams amgx->jacobi_relaxation_factor = 0.9;
609e6f8f311SMark Adams amgx->gs_symmetric = PETSC_FALSE;
610e6f8f311SMark Adams amgx->print_grid_stats = PETSC_FALSE;
611e6f8f311SMark Adams amgx->verbose = PETSC_FALSE;
612e6f8f311SMark Adams amgx->rsrc_init = false;
613e6f8f311SMark Adams amgx->solve_state_init = false;
614e6f8f311SMark Adams
615e6f8f311SMark Adams s_count++;
616e6f8f311SMark Adams
617e6f8f311SMark Adams PetscCallCUDA(cudaGetDevice(&amgx->devID));
618e6f8f311SMark Adams if (s_count == 1) {
619e6f8f311SMark Adams PetscCallAmgX(AMGX_initialize());
620e6f8f311SMark Adams PetscCallAmgX(AMGX_initialize_plugins());
621e6f8f311SMark Adams PetscCallAmgX(AMGX_register_print_callback(&print_callback));
622e6f8f311SMark Adams PetscCallAmgX(AMGX_install_signal_handler());
623e6f8f311SMark Adams }
624e6f8f311SMark Adams /* This communicator is not yet known to this system, so we duplicate it and make an internal communicator */
625e6f8f311SMark Adams PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &amgx->comm));
626e6f8f311SMark Adams PetscCallMPI(MPI_Comm_size(amgx->comm, &amgx->nranks));
627e6f8f311SMark Adams PetscCallMPI(MPI_Comm_rank(amgx->comm, &amgx->rank));
628e6f8f311SMark Adams
6293ba16761SJacob Faibussowitsch PetscCall(amgx_output_messages(amgx));
6303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
631e6f8f311SMark Adams }
632e6f8f311SMark Adams
633a22370e2Smarkadams4 /*@C
634e6f8f311SMark Adams PCAmgXGetResources - get AMGx's internal resource object
635e6f8f311SMark Adams
636cc4c1da9SBarry Smith Not Collective, No Fortran Support
637e6f8f311SMark Adams
638f1580f4eSBarry Smith Input Parameter:
639e6f8f311SMark Adams . pc - the PC
640e6f8f311SMark Adams
641e6f8f311SMark Adams Output Parameter:
642e6f8f311SMark Adams . rsrc_out - pointer to the AMGx resource object
643e6f8f311SMark Adams
644e6f8f311SMark Adams Level: advanced
645e6f8f311SMark Adams
646562efe2eSBarry Smith .seealso: [](ch_ksp), `PCAMGX`, `PC`, `PCGAMG`
647e6f8f311SMark Adams @*/
PCAmgXGetResources(PC pc,void * rsrc_out)648d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCAmgXGetResources(PC pc, void *rsrc_out)
649d71ae5a4SJacob Faibussowitsch {
650e6f8f311SMark Adams PC_AMGX *amgx = (PC_AMGX *)pc->data;
651e6f8f311SMark Adams
652e6f8f311SMark Adams PetscFunctionBegin;
653e6f8f311SMark Adams if (!amgx->rsrc_init) {
654e6f8f311SMark Adams // Read configuration file
655e6f8f311SMark Adams PetscCallAmgX(AMGX_config_create(&amgx->cfg, amgx->cfg_contents.c_str()));
656e6f8f311SMark Adams PetscCallAmgX(AMGX_resources_create(&amgx->rsrc, amgx->cfg, &amgx->comm, 1, &amgx->devID));
657e6f8f311SMark Adams amgx->rsrc_init = true;
658e6f8f311SMark Adams }
659e6f8f311SMark Adams *static_cast<AMGX_resources_handle *>(rsrc_out) = amgx->rsrc;
6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
661e6f8f311SMark Adams }
662