Page MenuHomec4science

micpsolver_base.hpp
No OneTemporary

File Metadata

Created
Wed, May 15, 20:17

micpsolver_base.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_MICPSOLVER_MICPSOLVERBASE_HPP
#define SPECMICP_MICPSOLVER_MICPSOLVERBASE_HPP
//! \file micpsolver_base.hpp Base class for micpsolver
#include <memory>
#include "common.hpp"
#include "micpsolver_structs.hpp"
#include "utils/log.hpp"
#include "utils/options_handler.hpp"
#include "utils/perfs_handler.hpp"
namespace specmicp {
namespace micpsolver {
//! \brief Base class for an MiCP solver
//!
//! Handle the program, options and performance
//!
template <class program_t>
class MiCPSolverBaseProgram:
public OptionsHandler<MiCPSolverOptions>,
public PerformanceHandler<MiCPPerformance>
{
public:
using program_ptr = std::shared_ptr<program_t>;
using OptionsHandler<MiCPSolverOptions>::get_options;
using PerformanceHandler<MiCPPerformance>::get_perfs;
MiCPSolverBaseProgram():
OptionsHandler<MiCPSolverOptions>(),
PerformanceHandler<MiCPPerformance>()
{}
//! \brief Constructor
//!
//! \param prog shared_ptr to the program
MiCPSolverBaseProgram(program_ptr prog):
OptionsHandler<MiCPSolverOptions>(),
PerformanceHandler<MiCPPerformance>(),
m_program(prog)
{}
//! \brief Return the program
//!
//! This is due to the templating, allow to do get_program()->whatever();
program_ptr& get_program() {return m_program;}
//! \brief Return the number of equations
index_t get_neq() const {return m_program->total_variables();}
//! \brief Return the number of equations corresponding to the free variables (size of G)
index_t get_neq_free() const {return m_program->nb_free_variables();}
//! \brief Compute the residuals, store it in r
//!
//! \param[in] x the variables
//! \param[out] r vector to store the residuals (must be of the same size as x)
void compute_residuals(const Vector& x, Vector& r)
{
m_program->get_residuals(x, r);
get_perfs().nb_call_residuals += 1;
}
//! \brief Compute the jacobian
//!
//! Assumes that the residual have been computed before
//!
//! \param[in] x the variables
//! \param[out] jacobian the jacobian
void compute_jacobian(Vector& x, Matrix& jacobian)
{
m_program->get_jacobian(x, jacobian);
get_perfs().nb_call_jacobian += 1;
}
//! \brief Compute the factors to scale the jacobian
//!
//! \param[in] jacobian the jacobian to scale (from the reformulated problem)
//! \param[in] residuals the residuals corresponding to the jacobian
//! \param[out] rscaler scaling factors of the rows
//! \param[out] cscaler scaling factors of the columns
void scaling_jacobian(
const Matrix& jacobian,
const Vector& residuals,
Vector& rscaler,
Vector& cscaler);
//! \brief Check for convergence
//!
//! \param nb_iterations the number of iterations
//! \param update the update taken at the previous iteration
//! \param solution the current solution
//! \param residuals current residuals
//! \param may_have_converged a boolean given by the user to indicates if the problem has converged yet
//! \return a MiCPSolverReturnCode describing the state of the algorithm
MiCPSolverReturnCode check_convergence(int nb_iterations,
const Vector& update,
const Vector& solution,
const Vector &residuals,
bool may_have_converged = true);
//! \brief Reformulate the jacobian using the cck function
//!
//! \param x the variables
//! \param r the residuals
//! \param jacobian the jacobian
void reformulate_jacobian_cck(
const Vector& x,
const Vector& r,
Matrix& jacobian
);
//! \brief Compute the norm of the update
//!
//! \tparam p order of the norm (1, 2, ..., Eigen::Infinity)
//! \param update the solution of the linear system
//! \param solution the solution (variables)
template <int p>
double norm_update(const Vector& update,
const Vector& solution) const {
return (update.array().abs()/(solution.array().max(1))
).matrix().template lpNorm<p>();
}
//! \brief Project variables on the feasible set
//!
//! \param x the variables
void projection(Vector& x);
//! \brief Return the step corrected step length if it is too long
//!
//! \param update the solution of the linear system
scalar_t is_step_too_long(Vector& update);
protected:
//! \brief Set the program
//!
//! \param program_ptr shared_ptr to the program
program_ptr& set_program(program_ptr theprogram) {
m_program = theprogram;
}
private:
program_ptr m_program;
};
} // end namespace micpsolver
} // end namespace specmicp
// ----------------------------------- //
// Implementation //
// ----------------------------------- //
#include "micpsolver_base.inl"
#endif // SPECMICP_MICPSOLVER_MICPSOLVERBASE_HPP

Event Timeline