ddemidov / amgcl

C++ library for solving large sparse linear systems with algebraic multigrid method
http://amgcl.readthedocs.org/
MIT License
728 stars 111 forks source link

Persisting warning: AMGCL WARNING: unknown parameter relaxation #90

Closed klausbu closed 3 years ago

klausbu commented 6 years ago

I tried to apply relaxation:

typedef amgcl::make_solver<
    amgcl::amg<
        Backend,
        amgcl::runtime::coarsening::wrapper,
        amgcl::runtime::relaxation::wrapper
        >,
    amgcl::runtime::solver::wrapper<Backend>
    > Solver;

    boost::property_tree::ptree prm;

    prm.put("solver.type", "lgmres");
    prm.put("solver.tol", tolerance_); 
    prm.put("solver.maxiter", 100);
    prm.put("precond.coarsening.type", "smoothed_aggregation");
    prm.put("precond.relaxation.type", "spai1"); 
    prm.put("precond.relaxation.damping", 0.7);

I replaced the last two lines also with

    prm.put("precond.relax.type", "spai1"); 
    prm.put("precond.relax.damping", 0.7);

as shows in the example in the paper "SUBDOMAIN DEFLATION COMBINED WITH LOCAL AMG: A CASE STUDY USING AMGCL LIBRARY"

But I get in every case a warning "AMGCL WARNING: unknown parameter relaxation"

What's the correct syntax?

ddemidov commented 6 years ago

The correct syntax is

    prm.put("precond.relax.type", "spai1"); 
    prm.put("precond.relax.damping", 0.7);

I am sorry for somewhat inconsistent naming scheme here. Here is how you can check the parameter structure:

https://github.com/ddemidov/amgcl/blob/5522c99408881a87ee1e6fadf8ce989cc3a66c0d/amgcl/make_solver.hpp#L67-L68

https://github.com/ddemidov/amgcl/blob/5522c99408881a87ee1e6fadf8ce989cc3a66c0d/amgcl/amg.hpp#L98-L146

etc etc. Each parent component includes its children's params into its own params, and when a runtime interface is used, the resulting structure is mapped to boost::property_tree during initialization.

The fact that you are still getting AMGCL WARNING: unknown parameter relaxation means that you probably still have precond.relaxation somewhere. The warning is issued when there is an element of boost::property_tree that the selected component does not know how to handle.

Also, just a note: spai1 usually does not pay off, since the setup is too expensive, and it does not result in any faster convergence. ilu0 may work better if spai0 convergence is not good.

ddemidov commented 6 years ago

Also, spai1 does not have damping parameter. In fact, its params struct is empty:

https://github.com/ddemidov/amgcl/blob/5522c99408881a87ee1e6fadf8ce989cc3a66c0d/amgcl/relaxation/spai1.hpp#L60

klausbu commented 6 years ago

Using AMGCL, what comes closest to a "Parallel Incomplete Cholesky" preconditioner?

The combination CG + IC is usually the clear winner, when it comes to solving the p equation of a CFD problem.

ddemidov commented 6 years ago

ilu0 would be the closest amgcl has to incomplete Cholesky. It is OpenMP-parallel (and CUDA backend uses ilu0 implementation from CUSPARSE). You can either use it directly, with

amgcl::make_solver<
  amgcl::relaxation::as_preconditioner<Backend, amgcl::relaxation::ilu0>,
  amgcl::solver::cg<Backend>
  >

or as a smoother within AMG:

amgcl::make_solver<
  amgcl::amg<
    Backend,
    amgcl::coarsening::smoothed_aggregation,
    amgcl::relaxation::ilu0
    >,
  amgcl::solver::cg<Backend>
  >
klausbu commented 6 years ago

How can I replace the amg with the pure ilu0 preconditioner in the mixed precision setup?

    // Preconditioner is in single precision:
    typedef
        amgcl::amg<
            amgcl::backend::builtin<float>,
            amgcl::coarsening::smoothed_aggregation,
            amgcl::relaxation::ilu0
            >
        Precond;

    // Solver is in double precision:
    typedef
        amgcl::solver::cg<
            amgcl::backend::builtin<double>
            >
        Solver;

    auto A_f = std::tie(n, amgcl_rows, amgcl_cols, amgcl_vals_fp32);
    auto A_d = std::tie(n, amgcl_rows, amgcl_cols, amgcl_vals);

    amgcl::profiler<> prof;

    prof.tic("setup");
    Solver S(n);
    Precond P(A_f);
    prof.toc("setup");

    int iters;
    double error;
    prof.tic("solve");
    std::tie(iters, error) = S(A_d, P, amgcl_rhs, amgcl_psi);
    prof.toc("solve");

Is it possible to move just the float operations to the GPU in a following step?

ddemidov commented 6 years ago

How can I replace the amg with the pure ilu0 preconditioner in the mixed precision setup?

You should just replace the definition of Precond to:

// Preconditioner is in single precision:
typedef
    amgcl::relaxation::as_preconditioner<
        amgcl:backend::builtin<float>,
        amgcl::relaxation::ilu0
        >
    Precond;

Is it possible to move just the float operations to the GPU in a following step?

This is not possible right now (you can not mix different backends together, only precision within same backend). Also, I don't think that would be a good idea, since it would involve heavy traffic between solver and preconditioner via slow PCI bus, which means it would be slower that either purely CPU or purely GPU approach.