Page MenuHomec4science

sparse_gmres.hpp
No OneTemporary

File Metadata

Created
Sun, May 5, 12:55

sparse_gmres.hpp

/*-------------------------------------------------------------------------------
Copyright (c) 2014,2015 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_SPARSEGMRES_HPP
#define SPECMICP_SPARSESOLVERS_SPARSEGMRES_HPP
#include "sparse_solver_structs.hpp"
#include "sparse_solver_base.hpp"
// needed for eigen 3.2.2
#include <iostream>
#include <Eigen/IterativeSolvers>
namespace specmicp {
namespace sparse_solvers {
template <class MatrixT, class DerivedR, class DerivedS>
class SparseSolverGMRES: public SparseSolverBase<MatrixT, DerivedR, DerivedS>
{
using SolverT = Eigen::GMRES<MatrixT, Eigen::IncompleteLUT<typename MatrixT::Scalar>> ;
public:
void analyse_pattern(MatrixT& jacobian)
{
}
//! \brief Decompose the jacboian
SparseSolverReturnCode decompose(MatrixT& jacobian)
{
m_solver.compute(jacobian);
if (m_solver.info() != Eigen::Success)
{
return SparseSolverReturnCode::FailedDecomposition;
}
return SparseSolverReturnCode::Success;
}
//! \brief Solve the problem
SparseSolverReturnCode solve(
const DerivedR& residuals,
DerivedS& solution
)
{
solution = m_solver.solve(-residuals);
if (m_solver.info() != Eigen::Success)
{
return SparseSolverReturnCode::FailedSystemSolving;
}
return SparseSolverReturnCode::Success;
}
SparseSolverReturnCode solve_scaling(
const DerivedR& residuals,
const DerivedR& scaling,
DerivedS& solution
)
{
solution = m_solver.solve(scaling.asDiagonal()*(-residuals));
if (m_solver.info() != Eigen::Success)
{
return SparseSolverReturnCode::FailedSystemSolving;
}
return SparseSolverReturnCode::Success;
}
private:
SolverT m_solver;
};
} // end namespace sparse_solvers
} // end namespace specmicp
#endif //SPECMICP_SPARSESOLVERS_SPARSEGMRES_HPP

Event Timeline