Page MenuHomec4science

sparse_solver.hpp
No OneTemporary

File Metadata

Created
Thu, Nov 7, 01:59

sparse_solver.hpp

/* =============================================================================
Copyright (c) 2014 - 2016
F. Georget <fabieng@princeton.edu> Princeton University
All rights reserved.
Redistribution and use in source and binary forms, with or without modification,
are permitted provided that the following conditions are met:
1. Redistributions of source code must retain the above copyright notice,
this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright notice,
this list of conditions and the following disclaimer in the documentation
and/or other materials provided with the distribution.
3. Neither the name of the copyright holder nor the names of its
contributors may be used to endorse or promote products derived from this
software without specific prior written permission.
* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. *
============================================================================= */
#ifndef SPECMICP_SPARSESOLVERS_SPARSESOLVERS_HPP
#define SPECMICP_SPARSESOLVERS_SPARSESOLVERS_HPP
/*!
\file sparse_solvers/sparse_solver.hpp
\brief Wrapper around Eigen sparse solvers
The wrapper allow to choose the solver at runtime
\code{.cpp}
auto solver = get_sparse_solver<MatrixT,Vector,Vector>(SparseSolver::SparseLU);
auto retcode = solver.decompose(Jacobian);
if (retcode != Success) {
throw std::runtime_error("Error during the decomposition");
}
retcode = solver.solver(residuals, solution);
if (retcode != Success) {
throw std::runtime_error("Error : can't solve the system");
}
\endcode
*/
#include <memory>
#include "specmicp_common/types.hpp"
#include "specmicp_common/compat.hpp"
#include "specmicp_common/eigen/incl_eigen_sparse_core.hpp"
// The main solvers available
#include "sparse_qr.hpp"
#include "sparse_lu.hpp"
#include "sparse_bicgstab.hpp"
// only include gmres if it is available
// from the Eigen library unsuported modules
#ifdef EIGEN_UNSUPPORTED_FOUND
#include "sparse_gmres.hpp"
#endif
namespace specmicp {
//! \namespace specmicp::sparse_solvers
//! \brief Wrappers around the Eigen sparse solver
namespace sparse_solvers {
//! \brief A pointer to a sparse solver
//!
//! \tparam MatrixT type of the matrix
//! \tparam DerivedR type of the residuals
//! \tparam DerivedS type of the solution
template <class MatrixT, class DerivedR, class DerivedS>
using SparseSolverPtr = std::unique_ptr<
SparseSolverBase<MatrixT, DerivedR, DerivedS>>;
//! \brief Return a sparse solver of type 'solver_type'
//!
//! \tparam MatrixT type of the matrix
//! \tparam DerivedR type of the residuals
//! \tparam DerivedS type of the solution
template <class MatrixT, class DerivedR, class DerivedS>
SparseSolverPtr<MatrixT, DerivedR, DerivedS>
get_sparse_solver(SparseSolver solver_type)
{
SparseSolverPtr<MatrixT, DerivedR, DerivedS> sparse_solver;
switch (solver_type) {
case SparseSolver::SparseQR:
sparse_solver = make_unique<SparseSolverQR<MatrixT, DerivedR, DerivedS>>();
break;
#ifdef EIGEN_UNSUPPORTED_FOUND
case SparseSolver::GMRES:
sparse_solver = make_unique<SparseSolverGMRES<MatrixT, DerivedR, DerivedS>>();
break;
#endif
case SparseSolver::SparseLU:
sparse_solver = make_unique<SparseSolverLU<MatrixT, DerivedR, DerivedS>>();
break;
case SparseSolver::BiCGSTAB:
sparse_solver = make_unique<SparseSolverBiCGSTAB<MatrixT, DerivedR, DerivedS>>();
break;
}
return sparse_solver;
}
} // end namespace sparse_solvers
} // end namespace specmicp
#endif //SPECMICP_SPARSESOLVERS_SPARSESOLVERS_HPP

Event Timeline