Page MenuHomec4science

simplex_solver.hpp
No OneTemporary

File Metadata

Created
Mon, May 13, 19:19

simplex_solver.hpp

/* =============================================================================
Copyright (c) 2014-2017 F. Georget <fabieng@princeton.edu> Princeton University
Copyright (c) 2017-2019 F. Georget <fabien.georget@epfl.ch> EPFL
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_COMMON_SIMPLEX_SOLVER_HPP
#define SPECMICP_COMMON_SIMPLEX_SOLVER_HPP
//! \file simplex/simplex_solver.hpp
//! \brief Declaration of the simplex solver
#include "specmicp_common/types.hpp"
#include <Eigen/LU>
namespace specmicp {
//! \namespace specmicp::simplex
//! \brief A simplex solver
namespace simplex {
//! \brief The return code of the simplex algorithm
enum class SimplexSolverReturnCode
{
Error = -10, //!< Generic error
Unbounded = -1, //!< Problem is unbounded
NotConverged = 0, //!< Solver has not converged yet
Success = 1 //!< Solver has converged
};
//! \brief The Simplex solver
//!
//! Using the Revised Simplex algorithm
//! Ref: Nocedal & Wright (2006)
template<class SimplexProgram>
class SimplexSolver {
using MatrixType = Matrix; // has to be for now
using BasisType = Eigen::Matrix<index_t, Eigen::Dynamic, 1>;//
public:
SimplexSolver(SimplexProgram& program):
m_program(program)
{}
//! \brief Solve the linear system
//!
//! Return SimplexSolverReturnCode::Success if converged
SimplexSolverReturnCode solve();
//! \brief Solve one step of the simplex method
//!
//! Return SimplexSolverReturnCode::Success if converged
SimplexSolverReturnCode one_step();
//! \brief Return the solution
//!
//! \param[out] val vector of the variables
//! \return the value of the objective function
scalar_t get_solution(Vector& val);
private:
//! \brief Initialize the data structures and fill values
void initialize();
//! \brief The pricing step, how to improve the minimization
SimplexSolverReturnCode pricing_step();
//! \brief Check convergence
SimplexSolverReturnCode check_objective();
//! \brief Select the variable to enter the basis
index_t select_q();
//! \brief Needed to select which varible to exit
void solve_for_d(index_t q);
//! \brief Check if the problem is unbounded
SimplexSolverReturnCode check_unboundness();
//! \brief Select the variable to exit the basis
std::pair<index_t, scalar_t> select_p();
//! \brief Update the variables vector
void update_xB(scalar_t xq);
//! \brief Modify the basis
void modify_basis(index_t p, index_t q, scalar_t xq);
//! \brief Start the simplex algorithm
void compute_start_step();
SimplexProgram& m_program; //!< Reference to the program
BasisType m_B; //!< Basis of non-zero variables
BasisType m_N; //!< Basis of zero variables
Matrix m_mB; //!< The non-zero part of the equalities
Vector m_xB; //!< The non-zero variables
Vector m_sN; //!< Gradient of the objective function
Vector m_cB; //!< Coefficient of non-zero term the objective functions
Vector m_db; //!< Column of the equalities to switch in
Vector m_d; //!< Vector used to find the varaible to set to zero
index_t m_previous_in_basis {-1};
Eigen::PartialPivLU<MatrixType> m_LU; //!< LU decomposition
};
}
}
#endif // SPECMICP_COMMON_SIMPLEX_SOLVER_HPP
#ifndef SPECMICP_COMMON_SIMPLEX_SOLVER_INL
#include "specmicp_common/simplex/simplex_solver.inl"
#endif

Event Timeline