Page MenuHomec4science

simplex_solver.inl
No OneTemporary

File Metadata

Created
Sat, Feb 22, 06:46

simplex_solver.inl

/* =============================================================================
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_INL
#define SPECMICP_COMMON_SIMPLEX_SOLVER_INL
// for syntax coloring
#ifndef SPECMICP_COMMON_SIMPLEX_SOLVER_HPP
#include "simplex_solver.hpp"
#endif //SPECMICP_COMMON_SIMPLEX_SOLVER_HPP
#include "specmicp_common/eigen/vector_checker.hpp"
#include "specmicp_common/log.hpp"
namespace specmicp {
namespace simplex {
template<class SimplexProgram>
void SimplexSolver<SimplexProgram>::initialize()
{
auto m = m_program.size_B();
auto p = m_program.size_N();
m_B.setZero(m);
for (auto i=0; i<m; ++i) {
m_B(i)=double(i);
}
m_N.setZero(p);
for (auto i=0; i<p; ++i) {
m_N(i)=double(m+i);
}
m_xB.setZero(m);
m_program.get_starting_guess(m_xB);
m_sN.setZero(p);
m_cB.setZero(m);
m_program.get_cB(m_cB);
m_d.setZero(m);
m_db.setZero(m);
m_mB.setIdentity(m,m);
}
template<class SimplexProgram>
SimplexSolverReturnCode SimplexSolver<SimplexProgram>::solve()
{
initialize();
SimplexSolverReturnCode ret = SimplexSolverReturnCode::NotConverged;
index_t cnt = 0;
while (ret == SimplexSolverReturnCode::NotConverged and cnt<50) {
ret = one_step();
++cnt;
}
return ret;
}
template<class SimplexProgram>
SimplexSolverReturnCode SimplexSolver<SimplexProgram>::one_step()
{
compute_start_step();
SimplexSolverReturnCode ret = pricing_step();
if (ret != SimplexSolverReturnCode::NotConverged) {
return ret;
}
ret = check_objective();
if (ret != SimplexSolverReturnCode::NotConverged) {
return ret;
}
index_t q = select_q();
if (q==-1) {
return SimplexSolverReturnCode::Error;
}
solve_for_d(q);
ret = check_unboundness();
if (ret != SimplexSolverReturnCode::NotConverged) {
return ret;
}
std::pair<index_t, scalar_t> pandxq = select_p();
if (pandxq.first == -1) {
return SimplexSolverReturnCode::Error;
}
update_xB(pandxq.second);
modify_basis(pandxq.first, q, pandxq.second);
return SimplexSolverReturnCode::NotConverged;
}
template<class SimplexProgram>
void SimplexSolver<SimplexProgram>::compute_start_step()
{
m_LU.compute(m_mB);
}
template<class SimplexProgram>
SimplexSolverReturnCode SimplexSolver<SimplexProgram>::pricing_step()
{
Vector lamb = m_LU.transpose().solve(m_cB);
if (vector_checker::is_finite(lamb)) {
m_program.compute_sN(m_N, lamb, m_sN);
return SimplexSolverReturnCode::NotConverged;
} else {
ERROR << "Error in Simplex solver: lambda is non finite";
return SimplexSolverReturnCode::Error;
}
}
template <class SimplexProgram>
SimplexSolverReturnCode SimplexSolver<SimplexProgram>::check_objective()
{
// bool not_converged_yet = false;
// for (index_t r=0;r<m_program.size_N();++r) {
// if (std::isfinite(m_sN(r)) and m_sN(r) < 0.0) {
// not_converged_yet = true;
// }
// }
// if (not_converged_yet) {
// return SimplexSolverReturnCode::NotConverged;
// } else {
// return SimplexSolverReturnCode::Success;
// }
if (vector_checker::all(m_sN) >= 0.0) {
return SimplexSolverReturnCode::Success;
}
return SimplexSolverReturnCode::NotConverged;
}
template <class SimplexProgram>
index_t SimplexSolver<SimplexProgram>::select_q()
{
index_t q=-1;
scalar_t current_min = +infinity();
for (index_t k=0; k<m_program.size_N();++k) {
if (std::isfinite(m_sN(k)) and m_sN(k)<current_min and m_N(k) != m_previous_in_basis and m_N(k)>m_program.size_B()) {
q = k;
current_min = m_sN(k);
}
}
return q;
}
template <class SimplexProgram>
void SimplexSolver<SimplexProgram>::solve_for_d(index_t q)
{
m_program.fill_Aq(m_N, m_db, q);
m_d = m_LU.solve(m_db);
}
template <class SimplexProgram>
SimplexSolverReturnCode SimplexSolver<SimplexProgram>::check_unboundness()
{
if (vector_checker::all(m_d) <= 0.0) {
return SimplexSolverReturnCode::Unbounded;
}
return SimplexSolverReturnCode::NotConverged;
}
template <class SimplexProgram>
std::pair<index_t, scalar_t> SimplexSolver<SimplexProgram>::select_p()
{
index_t p = -1;
scalar_t x_q = +infinity();
for (index_t i=0;i<m_program.size_B();++i) {
const scalar_t comp = m_xB(i)/m_d(i);
if (std::isfinite(comp) and comp < x_q) {
x_q = comp;
p = i;
}
}
return std::pair<index_t, scalar_t>(p, x_q);
}
template <class SimplexProgram>
void SimplexSolver<SimplexProgram>::update_xB(scalar_t xq)
{
m_xB -= m_d*xq;
}
template <class SimplexProgram>
void SimplexSolver<SimplexProgram>::modify_basis(index_t p, index_t q, scalar_t xq)
{
// B+
m_mB.col(p) = m_db;
//update cB
m_cB(p) = m_program.cN_q(m_N, q);
const auto tmp = m_B(p);
m_B(p)=m_N(q);
m_N(q) = tmp;
m_previous_in_basis=tmp;
m_xB(p) = xq;
}
template <class SimplexProgram>
scalar_t SimplexSolver<SimplexProgram>::get_solution(Vector& x)
{
x = Vector::Zero(m_program.nb_variables());
for (index_t i=0; i<m_program.size_B(); ++i) {
auto r = m_B(i);
x(r) = m_xB(i);
}
return m_cB.dot(m_xB);
}
}
}
#endif // SPECMICP_COMMON_SIMPLEX_SOLVER_INL

Event Timeline