mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Let amgcl use runtime parameters via JSON file
This commit is contained in:
parent
a23d881817
commit
3c1bfeb72f
@ -27,6 +27,9 @@
|
||||
#include <opm/simulators/linalg/bda/BdaResult.hpp>
|
||||
#include <opm/simulators/linalg/bda/amgclSolverBackend.hpp>
|
||||
|
||||
#include <boost/property_tree/json_parser.hpp>
|
||||
|
||||
|
||||
namespace bda
|
||||
{
|
||||
|
||||
@ -69,11 +72,12 @@ void amgclSolverBackend<block_size>::initialize(int N_, int nnz_, int dim, doubl
|
||||
cusparseCreate(&bprm.cusparse_handle);
|
||||
#endif
|
||||
|
||||
// set amgcl parameters
|
||||
prm.precond.damping = 0.9;
|
||||
prm.solver.maxiter = maxit;
|
||||
prm.solver.tol = tolerance;
|
||||
prm.solver.verbose = (verbosity >= 2);
|
||||
// read amgcl parameters via json file
|
||||
std::ifstream file("amgcl_options.json");
|
||||
boost::property_tree::read_json(file, prm);
|
||||
|
||||
// print json parameters
|
||||
boost::property_tree::write_json(std::cout, prm);
|
||||
|
||||
initialized = true;
|
||||
} // end initialize()
|
||||
@ -143,7 +147,12 @@ void amgclSolverBackend<block_size>::solve_system(double *b, WellContributions &
|
||||
|
||||
try {
|
||||
// create matrix object
|
||||
#if AMGCL_CUDA
|
||||
auto A = std::tie(N, A_rows, A_cols, A_vals);
|
||||
#else
|
||||
auto Atmp = std::tie(N, A_rows, A_cols, A_vals);
|
||||
auto A = amgcl::adapter::block_matrix<dmat_type>(Atmp);
|
||||
#endif
|
||||
|
||||
// create solver and construct preconditioner
|
||||
// don't reuse this unless the preconditioner can be reused
|
||||
@ -174,7 +183,7 @@ void amgclSolverBackend<block_size>::solve_system(double *b, WellContributions &
|
||||
auto X = amgcl::make_iterator_range(x_ptr, x_ptr + N / block_size);
|
||||
|
||||
// actually solve
|
||||
std::tie(iters, error) = solve(A, B, X);
|
||||
std::tie(iters, error) = solve(B, X);
|
||||
#endif
|
||||
|
||||
} catch (const std::exception& ex) {
|
||||
|
@ -26,6 +26,8 @@
|
||||
#include <opm/simulators/linalg/bda/BdaSolver.hpp>
|
||||
#include <opm/simulators/linalg/bda/WellContributions.hpp>
|
||||
|
||||
#include <boost/property_tree/ptree.hpp>
|
||||
|
||||
#include <amgcl/amg.hpp>
|
||||
#include <amgcl/backend/builtin.hpp>
|
||||
#include <amgcl/backend/cuda.hpp>
|
||||
@ -36,6 +38,8 @@
|
||||
#include <amgcl/relaxation/ilu0.hpp>
|
||||
#include <amgcl/solver/bicgstab.hpp>
|
||||
|
||||
#include <amgcl/preconditioner/runtime.hpp>
|
||||
|
||||
#include <amgcl/value_type/static_matrix.hpp>
|
||||
|
||||
#define AMGCL_CUDA 0
|
||||
@ -72,15 +76,7 @@ class amgclSolverBackend : public BdaSolver<block_size>
|
||||
typedef amgcl::backend::builtin<dmat_type> Backend;
|
||||
#endif
|
||||
|
||||
// choose linear solver components
|
||||
typedef amgcl::relaxation::as_preconditioner<Backend, amgcl::relaxation::ilu0> Precond;
|
||||
typedef amgcl::solver::bicgstab<Backend> IterSolver;
|
||||
|
||||
#if AMGCL_CUDA
|
||||
typedef amgcl::make_solver<Precond, IterSolver> Solver;
|
||||
#else
|
||||
typedef amgcl::make_block_solver<Precond, IterSolver> Solver;
|
||||
#endif
|
||||
typedef amgcl::make_solver<amgcl::runtime::preconditioner<Backend>, amgcl::runtime::solver::wrapper<Backend> > Solver;
|
||||
|
||||
private:
|
||||
// store matrix in CSR format
|
||||
@ -89,7 +85,7 @@ private:
|
||||
std::vector<double> x;
|
||||
std::once_flag print_info;
|
||||
|
||||
typename Solver::params prm;
|
||||
boost::property_tree::ptree prm;
|
||||
typename Backend::params bprm; // amgcl backend parameters, only used for cusparseHandle
|
||||
|
||||
/// Initialize GPU and allocate memory
|
||||
|
Loading…
Reference in New Issue
Block a user