Page MenuHomec4science

micpsolver.hpp
No OneTemporary

File Metadata

Created
Sun, Jun 30, 06:14

micpsolver.hpp

/*-------------------------------------------------------
- Module : micpsolver
- File : micpsolver.hpp
- Author : Fabien Georget
Copyright (c) 2014, Fabien 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:
* Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
* 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.
* Neither the name of the Princeton University 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 OWNER 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 SPECMIC_MICPSOLVER_MICPSOLVER_HPP
#define SPECMIC_MICPSOLVER_MICPSOLVER_HPP
#include "common.hpp"
#include "micpsolver_base.hpp"
#include "ncp_function.hpp"
#include "Eigen/QR"
//! \file micpsolver.hpp The MiCP solver
namespace specmicp {
namespace micpsolver {
//! \class MiCPSolver
//! \brief The MiCP Solver
//!
//! Solve
//! - \f$ G(u, v) = 0\f$
//! - \f$ 0 \leq v \perp H(u,v) \geq 0 \f$
//!
//! \tparam program_t a subclass of MiCPProg
//!
//! References :
//! - \cite Munson2001
//! - \cite Facchinei2003
//!
template <class program_t>
class MiCPSolver: public MiCPSolverBaseProgram<program_t>
{
public:
using base = MiCPSolverBaseProgram<program_t>;
using MiCPSolverBaseProgram<program_t>::get_options;
using MiCPSolverBaseProgram<program_t>::get_program;
using MiCPSolverBaseProgram<program_t>::get_perfs;
using MiCPSolverBaseProgram<program_t>::get_neq;
using MiCPSolverBaseProgram<program_t>::get_neq_free;
//! \brief Constructor
//!
//! \param prog smart pointer toward an instance of a program_t to solve
MiCPSolver(std::shared_ptr<program_t> prog):
MiCPSolverBaseProgram<program_t>(prog),
m_residuals(prog->total_variables()),
m_phi_residuals(Vector::Zero(prog->total_variables())),
m_jacobian(prog->total_variables(),prog ->total_variables())
{
}
// Merit function
// ##############
// Residuals and jacobian
// ----------------------
//! \brief Compute the residuals, use internal storage
//!
//! \param[in] x the variables
void compute_residuals(const Vector& x)
{
base::compute_residuals(x, m_residuals);
}
//! \brief Compute the jacobian
//!
//! Assumes that the residual have been computed before
//!
//! \param[in] x the variables
void compute_jacobian(Vector& x)
{
base::compute_jacobian(x, m_jacobian);
}
//! \brief Reformulation of the residuals
//!
//! Reformulate the problem - assumes that r, the residual, has been computed before
//!
//! \param[in] x the variables
//! \param[in] r the residuals
//! \param[out] r_phi a vector of size neq, which will contain the reformulated residuals
void reformulate_residuals(const Vector& x,
const Vector& r,
Vector& r_phi);
//! \brief Reformulation of the residuals *inplace*
//!
//! \param[in] x the variables
//! \param[in,out] r the residual, will contain the reformulated residuals
void reformulate_residuals_inplace(const Vector &x,
Vector &r);
//! \brief Reformulation of the jacobian
//!
//! r is the original vector of residuals (not reformulated)
void reformulate_jacobian(const Vector& x);
// Algorithm
// #########
//! \brief Solver the program using x as initial guess
//!
//! \param[in,out] x the initial guess, as output contains the solution (from the last iteration)
MiCPSolverReturnCode solve(Vector& x);
//! \brief Setup the residuals
//!
//! \param[in] x the current solution
void setup_residuals(const Eigen::VectorXd& x) {
compute_residuals(x);
reformulate_residuals(x, m_residuals, m_phi_residuals);
}
//! \brief Setup the jacobian
//!
//! \param[in] x the current solution
void setup_jacobian(Eigen::VectorXd& x) {
compute_jacobian(x);
reformulate_jacobian(x);
m_grad_phi = m_jacobian.transpose() * m_phi_residuals;
}
//! \brief Solve the Newton system
//!
//! The system will be scaled before being solved, may be expensive and not necessary.
//! Do disable the scaling, disable the corresponding options in #MiCPSolverOptions.
//! It assumes that the Newton system has been formed
//!
//! \param[out] update the update to apply to the solution
//!
//! \sa search_direction_calculation_no_scaling.
MiCPSolverReturnCode search_direction_calculation(Vector& update);
//! \brief Solve the Newton system - does not scale the jacobian
//!
//! The system will not be scaled.
//! It assumes that the Newton system has been formed
//!
//! \param[out] update the update to apply to the solution
//!
//! \sa search_direction_scaling
MiCPSolverReturnCode search_direction_calculation_no_scaling(Vector &update);
//! \brief Linesearch
//!
//! It is a backtracking linesearch.
//! If the correct option is set, this is a nonmonotone linesearch.
//!
//! \param[in] x the current solution
//!
//! References :
//! - \cite Dennis1983
//! - \cite Munson2001
//!
//! \sa #MiCPSolverOptions
int linesearch(Eigen::VectorXd& update, Eigen::VectorXd& x);
//! \brief Crashing
//!
//! This function improves, if possible, the initial guess
//!
//! \param[in] x the current solution
//!
//! Reference :
//! - \cite Munson2001
void crashing(Vector &x);
private:
// Residuals and jacobian
Eigen::VectorXd m_residuals; //!< The residuals
Eigen::VectorXd m_phi_residuals; //!< The reformulated residuals
Eigen::VectorXd m_grad_phi; //!< The gradient of the reformulated residuals
Eigen::MatrixXd m_jacobian; //!< The jacobian
std::vector<scalar_t> m_max_merits; //!< Contains the m best value of the merit function
bool m_gradient_step_taken; //!< True if the update was computed using the gradient
};
} // end namespace micpsolver
} // end namespace specmicp
// ###############//
// Implementation //
// ###############//
#include "micpsolver.inl"
#endif // SPECMIC_MICPSOLVER_MICPSOLVER_HPP

Event Timeline