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