Advanced Multi-Physics (AMP)
On-Line Documentation
BoomerAMGSolver.h
Go to the documentation of this file.
1#ifndef included_AMP_BoomerAMGSolver
2#define included_AMP_BoomerAMGSolver
3
4
5#include "AMP/matrices/Matrix.h"
6#include "AMP/solvers/SolverStrategy.h"
7#include "AMP/solvers/SolverStrategyParameters.h"
8#include "AMP/solvers/hypre/HypreSolver.h"
9
10namespace AMP::LinearAlgebra {
11class HypreMatrixAdaptor;
12}
13
14namespace AMP::Solver {
15
16
18
19
28{
29
30public:
35
42 explicit BoomerAMGSolver( std::shared_ptr<BoomerAMGSolverParameters> parameters );
43
48
49 std::string type() const override { return "BoomerAMGSolver"; }
50
52 static std::unique_ptr<SolverStrategy>
53 createSolver( std::shared_ptr<SolverStrategyParameters> solverStrategyParameters )
54 {
55 return std::make_unique<BoomerAMGSolver>( solverStrategyParameters );
56 }
57
63 void apply( std::shared_ptr<const AMP::LinearAlgebra::Vector> f,
64 std::shared_ptr<AMP::LinearAlgebra::Vector> u ) override;
65
75 void initialize( std::shared_ptr<const SolverStrategyParameters> parameters ) override;
76
85 void reset( std::shared_ptr<SolverStrategyParameters> params ) override;
86
87 void getFromInput( std::shared_ptr<const AMP::Database> db );
88
89private:
90 void setupBoomerAMG( void );
92
97 int d_max_levels = 10;
101 int d_num_paths = 1;
103 int d_nodal = 0;
118 int d_relax_type = 16;
129 int d_rap2 = 0;
131
132 double d_strong_threshold = 0.25;
133 double d_max_row_sum = 0.9;
134 double d_non_galerkin_tol = 0.0;
135 double d_trunc_factor = 0.0;
136 double d_agg_trunc_factor = 0.0;
139 double d_relax_weight = 0.0;
140 double d_outer_weight = 1.0;
142 double d_schwarz_weight = 0.0;
143
144 bool d_bSetupSolver = true;
145};
146} // namespace AMP::Solver
147
148#endif
void reset(std::shared_ptr< SolverStrategyParameters > params) override
std::string type() const override
Return the name of the solver.
static std::unique_ptr< SolverStrategy > createSolver(std::shared_ptr< SolverStrategyParameters > solverStrategyParameters)
static create routine that is used by SolverFactory
void getFromInput(std::shared_ptr< const AMP::Database > db)
void initialize(std::shared_ptr< const SolverStrategyParameters > parameters) override
BoomerAMGSolver(std::shared_ptr< BoomerAMGSolverParameters > parameters)
void apply(std::shared_ptr< const AMP::LinearAlgebra::Vector > f, std::shared_ptr< AMP::LinearAlgebra::Vector > u) override
std::shared_ptr< AMP::Solver::SolverStrategy > shared_ptr



Advanced Multi-Physics (AMP)
Oak Ridge National Laboratory
Idaho National Laboratory
Los Alamos National Laboratory
This page automatically produced from the
source code by doxygen
Last updated: Tue Mar 10 2026 13:06:41.
Comments on this page