Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F122191257
sparse_solver.hpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Wed, Jul 16, 12:06
Size
2 KB
Mime Type
text/x-c++
Expires
Fri, Jul 18, 12:06 (2 d)
Engine
blob
Format
Raw Data
Handle
27447364
Attached To
rSPECMICP SpecMiCP / ReactMiCP
sparse_solver.hpp
View Options
#ifndef SPECMICP_DFPMSOLVER_SPARSESOLVER_HPP
#define SPECMICP_DFPMSOLVER_SPARSESOLVER_HPP
#include <Eigen/Sparse>
#include "sparse_solver_structs.hpp"
namespace specmicp {
namespace dfpmsolver {
// header
// ======
//! \brief Solve a sparse linear system
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_sparse_linear_system(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update, SparseSolver solver);
//! \brief Solve a sparse linear system using the QR decomposition
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_sparseqr(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update);
//! \brief Solve a sparse linear system using the LU decomposition
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_sparselu(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update);
//! \brief Solve a sparse linear system using the BiCGSTAB algorithm
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_bicgstab(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update);
// Implementation
// ==============
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_sparse_linear_system(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update, SparseSolver solver)
{
if (solver == SparseSolver::SparseQR) {
return solve_linear_system_sparseqr(residual, jacobian, update);
} else if (solver == SparseSolver::SparseLU) {
return solve_linear_system_sparselu(residual, jacobian, update);
}
else {
return solve_linear_system_bicgstab(residual, jacobian, update);
}
}
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_sparseqr(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update)
{
Eigen::SparseQR<DerivedM, Eigen::COLAMDOrdering<int>> solver(jacobian);
if (solver.info() != Eigen::Success)
{
return -2;
}
update = solver.solve(-residual);
if (solver.info() != Eigen::Success)
{
return -1;
}
return 1;
}
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_sparselu(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update)
{
Eigen::SparseLU<DerivedM> solver(jacobian);
if (solver.info() != Eigen::Success)
{
return -2;
}
update = solver.solve(-residual);
if (solver.info() != Eigen::Success)
{
return -1;
}
return 1;
}
template <class DerivedV1, class DerivedV2, class DerivedM>
int solve_linear_system_bicgstab(DerivedV1& residual, DerivedM& jacobian, DerivedV2& update)
{
Eigen::BiCGSTAB<DerivedM, Eigen::IncompleteLUT<typename DerivedM::Scalar>> solver(jacobian);
update = solver.solve(-residual);
if (solver.info() != Eigen::Success)
{
return -1;
}
return 1;
}
} // end namespace dfpmsolver
} // end namespace specmicp
#endif //SPECMICP_DFPMSOLVER_SPARSESOLVER_HPP
Event Timeline
Log In to Comment