Lines Matching +full:- +full:std

7      pcimpl.h - private include file intended for use by all preconditioners
59 static const std::map<std::string, AmgXAMGMethod> AMGMethods;
60 static const std::map<std::string, AmgXSmoother> Smoothers;
61 static const std::map<std::string, AmgXSelector> Selectors;
62 static const std::map<std::string, AmgXCoarseSolver> CoarseSolvers;
63 static const std::map<std::string, AmgXAMGCycle> AMGCycles;
66 const std::map<std::string, AmgXAMGMethod> AmgXControlMap::AMGMethods = {
71 const std::map<std::string, AmgXSmoother> AmgXControlMap::Smoothers = {
87 const std::map<std::string, AmgXSelector> AmgXControlMap::Selectors = {
96 const std::map<std::string, AmgXCoarseSolver> AmgXControlMap::CoarseSolvers = {
101 const std::map<std::string, AmgXAMGCycle> AmgXControlMap::AMGCycles = {
130 std::string cfg_contents;
132 // Cached state for re-setup
164 static std::string amgx_output{};
177 if (amgx->verbose && !amgx_output.empty()) {
179 PetscCall(PetscPrintf(amgx->comm, "AMGX: %s", amgx_output.c_str()));
198 SETERRQ(amgx->comm, PETSC_ERR_LIB, "%s", msg); \
203 PCSetUp_AMGX - Prepares for the use of the AmgX preconditioner
207 . pc - the preconditioner context
217 PC_AMGX *amgx = (PC_AMGX *)pc->data;
218 Mat Pmat = pc->pmat;
225 // Non-sequential/MPI matrices must be adapted to extract the local matrix
226 bool partial_setup_allowed = (pc->setupcalled && pc->flag != DIFFERENT_NONZERO_PATTERN);
227 if (amgx->nranks > 1) {
229 PetscCall(MatMPIAIJGetLocalMat(Pmat, MAT_REUSE_MATRIX, &amgx->localA));
231 PetscCall(MatMPIAIJGetLocalMat(Pmat, MAT_INITIAL_MATRIX, &amgx->localA));
234 if (is_dev_ptrs) PetscCall(MatConvert(amgx->localA, MATSEQAIJCUSPARSE, MAT_INPLACE_MATRIX, &amgx->localA));
236 amgx->localA = Pmat;
240 PetscCall(MatSeqAIJCUSPARSEGetArrayRead(amgx->localA, &amgx->values));
242 PetscCall(MatSeqAIJGetArrayRead(amgx->localA, &amgx->values));
247 if (!amgx->rsrc_init) {
249 PetscCallAmgX(AMGX_config_create(&amgx->cfg, amgx->cfg_contents.c_str()));
250 PetscCallAmgX(AMGX_resources_create(&amgx->rsrc, amgx->cfg, &amgx->comm, 1, &amgx->devID));
251 amgx->rsrc_init = true;
254 PetscCheck(!amgx->solve_state_init, amgx->comm, PETSC_ERR_PLIB, "AmgX solve state initialisation already called.");
255 PetscCallAmgX(AMGX_matrix_create(&amgx->A, amgx->rsrc, AMGX_mode_dDDI));
256 PetscCallAmgX(AMGX_vector_create(&amgx->sol, amgx->rsrc, AMGX_mode_dDDI));
257 PetscCallAmgX(AMGX_vector_create(&amgx->rhs, amgx->rsrc, AMGX_mode_dDDI));
258 PetscCallAmgX(AMGX_solver_create(&amgx->solver, amgx->rsrc, AMGX_mode_dDDI, amgx->cfg));
259 amgx->solve_state_init = true;
265 PetscCall(MatGetRowIJ(amgx->localA, 0, PETSC_FALSE, PETSC_FALSE, &amgx->nLocalRows, &rowOffsets, &colIndices, &done));
266 PetscCheck(done, amgx->comm, PETSC_ERR_PLIB, "MatGetRowIJ was not successful");
267 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);
270 PetscCallCUDA(cudaMemcpy(&amgx->nnz, &rowOffsets[amgx->nLocalRows], sizeof(int), cudaMemcpyDefault));
272 amgx->nnz = rowOffsets[amgx->nLocalRows];
275 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);
278 std::vector<PetscInt> partitionOffsets(amgx->nranks + 1);
282 PetscCallMPI(MPI_Allgather(&amgx->nLocalRows, 1, MPIU_INT, partitionOffsets.data() + 1, 1, MPIU_INT, amgx->comm));
283 std::partial_sum(partitionOffsets.begin(), partitionOffsets.end(), partitionOffsets.begin());
286 amgx->nGlobalRows = partitionOffsets[amgx->nranks];
288 PetscCall(MatGetBlockSize(Pmat, &amgx->bSize));
290 // XXX Currently constrained to 32-bit indices, to be changed in the future
293 PetscCallAmgX(AMGX_distribution_create(&dist, amgx->cfg));
296 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));
297 PetscCallAmgX(AMGX_solver_setup(amgx->solver, amgx->A));
298 PetscCallAmgX(AMGX_vector_bind(amgx->sol, amgx->A));
299 PetscCallAmgX(AMGX_vector_bind(amgx->rhs, amgx->A));
302 PetscCall(MatRestoreRowIJ(amgx->localA, 0, PETSC_FALSE, PETSC_FALSE, &nlr, &rowOffsets, &colIndices, &done));
305 PetscCallAmgX(AMGX_matrix_replace_coefficients(amgx->A, amgx->nLocalRows, amgx->nnz, amgx->values, NULL));
306 PetscCallAmgX(AMGX_solver_resetup(amgx->solver, amgx->A));
310 PetscCall(MatSeqAIJCUSPARSERestoreArrayRead(amgx->localA, &amgx->values));
312 PetscCall(MatSeqAIJRestoreArrayRead(amgx->localA, &amgx->values));
319 PCApply_AMGX - Applies the AmgX preconditioner to a vector.
322 . pc - the preconditioner context
323 . b - rhs vector
326 . x - solution vector
332 PC_AMGX *amgx = (PC_AMGX *)pc->data;
348 PetscCallAmgX(AMGX_vector_upload(amgx->sol, amgx->nLocalRows, 1, x_));
349 PetscCallAmgX(AMGX_vector_upload(amgx->rhs, amgx->nLocalRows, 1, b_));
350 PetscCallAmgX(AMGX_solver_solve_with_0_initial_guess(amgx->solver, amgx->rhs, amgx->sol));
353 PetscCallAmgX(AMGX_solver_get_status(amgx->solver, &status));
355 PetscCheck(status != AMGX_SOLVE_FAILED, amgx->comm, PETSC_ERR_CONV_FAILED, "AmgX solver failed to solve the system! The error code is %d.", status);
356 PetscCallAmgX(AMGX_vector_download(amgx->sol, x_));
371 PC_AMGX *amgx = (PC_AMGX *)pc->data;
374 if (amgx->solve_state_init) {
375 PetscCallAmgX(AMGX_solver_destroy(amgx->solver));
376 PetscCallAmgX(AMGX_matrix_destroy(amgx->A));
377 PetscCallAmgX(AMGX_vector_destroy(amgx->sol));
378 PetscCallAmgX(AMGX_vector_destroy(amgx->rhs));
379 if (amgx->nranks > 1) PetscCall(MatDestroy(&amgx->localA));
381 amgx->solve_state_init = false;
387 PCDestroy_AMGX - Destroys the private context for the AmgX preconditioner
391 . pc - the preconditioner context
397 PC_AMGX *amgx = (PC_AMGX *)pc->data;
403 PetscCheck(amgx->rsrc != nullptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "s_rsrc == NULL");
404 PetscCallAmgX(AMGX_resources_destroy(amgx->rsrc));
406 PetscCallAmgX(AMGX_config_destroy(amgx->cfg));
409 PetscCallMPI(MPI_Comm_free(&amgx->comm));
411 PetscCallAmgX(AMGX_config_destroy(amgx->cfg));
413 s_count -= 1;
419 std::string map_reverse_lookup(const std::map<std::string, T> &map, const T &key)
429 PC_AMGX *amgx = (PC_AMGX *)pc->data;
435 amgx->cfg_contents = "config_version=2,";
436 amgx->cfg_contents += "determinism_flag=1,";
439 PetscCall(PetscOptionsBool("-pc_amgx_exact_coarse_solve", "AmgX AMG Exact Coarse Solve", "", amgx->exact_coarse_solve, &amgx->exact_coarse_solve, NULL));
440 if (amgx->exact_coarse_solve) amgx->cfg_contents += "exact_coarse_solve=1,";
442 amgx->cfg_contents += "solver(amg)=AMG,";
445 std::string def_amg_method = map_reverse_lookup(AmgXControlMap::AMGMethods, amgx->amg_method);
447 PetscCall(PetscOptionsString("-pc_amgx_amg_method", "AmgX AMG Method", "", option, option, MAX_PARAM_LEN, NULL));
449 amgx->amg_method = AmgXControlMap::AMGMethods.at(option);
450 amgx->cfg_contents += "amg:algorithm=" + std::string(option) + ",";
453 std::string def_amg_cycle = map_reverse_lookup(AmgXControlMap::AMGCycles, amgx->amg_cycle);
455 PetscCall(PetscOptionsString("-pc_amgx_amg_cycle", "AmgX AMG Cycle", "", option, option, MAX_PARAM_LEN, NULL));
457 amgx->amg_cycle = AmgXControlMap::AMGCycles.at(option);
458 amgx->cfg_contents += "amg:cycle=" + std::string(option) + ",";
461 std::string def_smoother = map_reverse_lookup(AmgXControlMap::Smoothers, amgx->smoother);
463 PetscCall(PetscOptionsString("-pc_amgx_smoother", "AmgX Smoother", "", option, option, MAX_PARAM_LEN, NULL));
465 amgx->smoother = AmgXControlMap::Smoothers.at(option);
466 amgx->cfg_contents += "amg:smoother(smooth)=" + std::string(option) + ",";
468 if (amgx->smoother == AmgXSmoother::JacobiL1 || amgx->smoother == AmgXSmoother::BlockJacobi) {
469 PetscCall(PetscOptionsScalar("-pc_amgx_jacobi_relaxation_factor", "AmgX AMG Jacobi Relaxation Factor", "", amgx->jacobi_relaxation_factor, &amgx->jacobi_relaxation_factor, NULL));
470 amgx->cfg_contents += "smooth:relaxation_factor=" + std::to_string(amgx->jacobi_relaxation_factor) + ",";
471 } else if (amgx->smoother == AmgXSmoother::GS || amgx->smoother == AmgXSmoother::MulticolorGS) {
472 PetscCall(PetscOptionsScalar("-pc_amgx_gs_symmetric", "AmgX AMG Gauss Seidel Symmetric", "", amgx->gs_symmetric, &amgx->gs_symmetric, NULL));
473 amgx->cfg_contents += "smooth:symmetric_GS=" + std::to_string(amgx->gs_symmetric) + ",";
477 std::string def_selector = map_reverse_lookup(AmgXControlMap::Selectors, amgx->selector);
479 PetscCall(PetscOptionsString("-pc_amgx_selector", "AmgX Selector", "", option, option, MAX_PARAM_LEN, NULL));
483 if (amgx->amg_method == AmgXAMGMethod::Classical) {
484 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);
485 amgx->cfg_contents += "amg:interpolator=D2,";
486 } else if (amgx->amg_method == AmgXAMGMethod::Aggregation) {
487 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");
489 amgx->selector = AmgXControlMap::Selectors.at(option);
490 amgx->cfg_contents += "amg:selector=" + std::string(option) + ",";
493 PetscCall(PetscOptionsInt("-pc_amgx_presweeps", "AmgX AMG Presweep Count", "", amgx->presweeps, &amgx->presweeps, NULL));
494 amgx->cfg_contents += "amg:presweeps=" + std::to_string(amgx->presweeps) + ",";
497 PetscCall(PetscOptionsInt("-pc_amgx_postsweeps", "AmgX AMG Postsweep Count", "", amgx->postsweeps, &amgx->postsweeps, NULL));
498 amgx->cfg_contents += "amg:postsweeps=" + std::to_string(amgx->postsweeps) + ",";
501 PetscCall(PetscOptionsInt("-pc_amgx_max_levels", "AmgX AMG Max Level Count", "", amgx->max_levels, &amgx->max_levels, NULL));
502 amgx->cfg_contents += "amg:max_levels=" + std::to_string(amgx->max_levels) + ",";
505 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));
506 amgx->cfg_contents += "amg:dense_lu_num_rows=" + std::to_string(amgx->dense_lu_num_rows) + ",";
509 PetscCall(PetscOptionsScalar("-pc_amgx_strength_threshold", "AmgX AMG Strength Threshold", "", amgx->strength_threshold, &amgx->strength_threshold, NULL));
510 amgx->cfg_contents += "amg:strength_threshold=" + std::to_string(amgx->strength_threshold) + ",";
513 PetscCall(PetscOptionsInt("-pc_amgx_aggressive_levels", "AmgX AMG Presweep Count", "", amgx->aggressive_levels, &amgx->aggressive_levels, NULL));
514 if (amgx->aggressive_levels > 0) amgx->cfg_contents += "amg:aggressive_levels=" + std::to_string(amgx->aggressive_levels) + ",";
517 std::string def_coarse_solver = map_reverse_lookup(AmgXControlMap::CoarseSolvers, amgx->coarse_solver);
519 PetscCall(PetscOptionsString("-pc_amgx_coarse_solver", "AmgX CoarseSolver", "", option, option, MAX_PARAM_LEN, NULL));
521 amgx->coarse_solver = AmgXControlMap::CoarseSolvers.at(option);
522 amgx->cfg_contents += "amg:coarse_solver=" + std::string(option) + ",";
525 amgx->cfg_contents += "amg:max_iters=1,";
528 PetscCall(PetscOptionsBool("-pc_amgx_print_grid_stats", "AmgX Print Grid Stats", "", amgx->print_grid_stats, &amgx->print_grid_stats, NULL));
530 if (amgx->print_grid_stats) amgx->cfg_contents += "amg:print_grid_stats=1,";
531 amgx->cfg_contents += "amg:monitor_residual=0";
534 PetscCall(PetscOptionsBool("-pc_amgx_verbose", "Enable output from AmgX", "", amgx->verbose, &amgx->verbose, NULL));
541 PC_AMGX *amgx = (PC_AMGX *)pc->data;
547 std::string output_cfg(amgx->cfg_contents);
548 std::replace(output_cfg.begin(), output_cfg.end(), ',', '\n');
555 PCAMGX - Interface to NVIDIA's AmgX algebraic multigrid
558 + -pc_amgx_amg_method <CLASSICAL,AGGREGATION> - set the AMG algorithm to use
559 . -pc_amgx_amg_cycle <V,W,F,CG> - set the AMG cycle type
560 . -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
561 . -pc_amgx_jacobi_relaxation_factor - set the relaxation factor for Jacobi smoothing
562 . -pc_amgx_gs_symmetric - enforce symmetric Gauss-Seidel smoothing (only applies if GS smoothing is selected)
563 . -pc_amgx_selector <SIZE_2,SIZE_4,SIZE_8,MULTI_PAIRWISE,PMIS,HMIS> - set the AMG coarse selector
564 . -pc_amgx_presweeps - set the number of AMG pre-sweeps
565 . -pc_amgx_postsweeps - set the number of AMG post-sweeps
566 . -pc_amgx_max_levels - set the maximum number of levels in the AMG level hierarchy
567 . -pc_amgx_strength_threshold - set the strength threshold for the AMG coarsening
568 . -pc_amgx_aggressive_levels - set the number of levels (from the finest) that should apply aggressive coarsening
569 . -pc_amgx_coarse_solver <DENSE_LU_SOLVER,NOSOLVER> - set the coarse solve
570 . -pc_amgx_print_grid_stats - output the AMG grid hierarchy to stdout
571 - -pc_amgx_verbose - enable AmgX output
587 pc->ops->apply = PCApply_AMGX;
588 pc->ops->setfromoptions = PCSetFromOptions_AMGX;
589 pc->ops->setup = PCSetUp_AMGX;
590 pc->ops->view = PCView_AMGX;
591 pc->ops->destroy = PCDestroy_AMGX;
592 pc->ops->reset = PCReset_AMGX;
593 pc->data = (void *)amgx;
596 amgx->selector = AmgXSelector::PMIS;
597 amgx->smoother = AmgXSmoother::BlockJacobi;
598 amgx->amg_method = AmgXAMGMethod::Classical;
599 amgx->coarse_solver = AmgXCoarseSolver::DenseLU;
600 amgx->amg_cycle = AmgXAMGCycle::V;
601 amgx->exact_coarse_solve = PETSC_TRUE;
602 amgx->presweeps = 1;
603 amgx->postsweeps = 1;
604 amgx->max_levels = 100;
605 amgx->strength_threshold = 0.5;
606 amgx->aggressive_levels = 0;
607 amgx->dense_lu_num_rows = 1;
608 amgx->jacobi_relaxation_factor = 0.9;
609 amgx->gs_symmetric = PETSC_FALSE;
610 amgx->print_grid_stats = PETSC_FALSE;
611 amgx->verbose = PETSC_FALSE;
612 amgx->rsrc_init = false;
613 amgx->solve_state_init = false;
617 PetscCallCUDA(cudaGetDevice(&amgx->devID));
625 PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &amgx->comm));
626 PetscCallMPI(MPI_Comm_size(amgx->comm, &amgx->nranks));
627 PetscCallMPI(MPI_Comm_rank(amgx->comm, &amgx->rank));
634 PCAmgXGetResources - get AMGx's internal resource object
639 . pc - the PC
642 . rsrc_out - pointer to the AMGx resource object
650 PC_AMGX *amgx = (PC_AMGX *)pc->data;
653 if (!amgx->rsrc_init) {
655 PetscCallAmgX(AMGX_config_create(&amgx->cfg, amgx->cfg_contents.c_str()));
656 PetscCallAmgX(AMGX_resources_create(&amgx->rsrc, amgx->cfg, &amgx->comm, 1, &amgx->devID));
657 amgx->rsrc_init = true;
659 *static_cast<AMGX_resources_handle *>(rsrc_out) = amgx->rsrc;