Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F60495759
micpsolver_base.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
Tue, Apr 30, 15:50
Size
5 KB
Mime Type
text/x-c++
Expires
Thu, May 2, 15:50 (2 d)
Engine
blob
Format
Raw Data
Handle
17361317
Attached To
rSPECMICP SpecMiCP / ReactMiCP
micpsolver_base.hpp
View Options
/*-------------------------------------------------------
- Module : micpsolver
- File : micpsolver_base.hpp
- Author : Fabien Georget
Copyright (c) 2014, Fabien Georget, Princeton University
---------------------------------------------------------*/
#ifndef SPECMICP_MICPSOLVER_MICPSOLVERBASE_HPP
#define SPECMICP_MICPSOLVER_MICPSOLVERBASE_HPP
//! \file micpsolver_base.hpp Base class for micpsolver
#include <memory>
#include "micpsolver_structs.hpp"
#include <Eigen/Core>
#include "utils/log.hpp"
namespace specmicp {
namespace micpsolver {
//! \brief Base class of an MiCP solver
//!
//! Handle options and performance structures.
class MiCPSolverBase
{
public:
//! \brief Return a const reference to the options used by the algorithm
const MiCPSolverOptions& get_options() const {return m_options;}
//! \brief Return a reference to the options used by the algorithm
MiCPSolverOptions& get_options() {return m_options;}
//! \brief Set the options
void set_options(const MiCPSolverOptions& options) {m_options = options;}
//! \brief Return a const reference to an instance of MiCPPerformance
const MiCPPerformance& get_performance() {return m_perf;}
protected:
//! \brief Return a reference to an instance of MiCPPerformance
MiCPPerformance& get_perf() {return m_perf;}
private:
MiCPSolverOptions m_options; //!< The options
MiCPPerformance m_perf; //!< The performance
};
//! \brief SubBase class of an MiCP solver, handle the program
//!
//!
template <class program_t>
class MiCPSolverBaseProgram: public MiCPSolverBase
{
public:
using program_ptr = std::shared_ptr<program_t>;
MiCPSolverBaseProgram():
MiCPSolverBase()
{}
MiCPSolverBaseProgram(program_ptr prog):
MiCPSolverBase(),
m_program(prog)
{}
//! \brief Return the program
program_ptr& get_program() {return m_program;}
//! \brief Return the number of equations
int get_neq() const {return m_program->total_variables();}
//! \brief Return the number of equations corresponding to the free variables (size of G)
int 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 Eigen::VectorXd& x, Eigen::VectorXd& r)
{
m_program->get_residuals(x, r);
get_perf().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(Eigen::VectorXd& x, Eigen::MatrixXd& jacobian)
{
m_program->get_jacobian(x, jacobian);
get_perf().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 Eigen::MatrixXd& jacobian,
const Eigen::VectorXd& residuals,
Eigen::VectorXd& rscaler,
Eigen::VectorXd& 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
//! @return a MiCPSolverReturnCode describing the state of the algorithm
MiCPSolverReturnCode check_convergence(int nb_iterations,
const Eigen::VectorXd &update,
const Eigen::VectorXd &solution,
const Eigen::VectorXd &residuals);
//! \brief Reformulate the jacobian using the cck function
void reformulate_jacobian_cck(
const Eigen::VectorXd& x,
const Eigen::VectorXd& r,
Eigen::MatrixXd& jacobian
);
//! \brief Compute the norm of the update
template <int p>
double norm_update(const Eigen::VectorXd& update,
const Eigen::VectorXd& solution) const {
return (update.array().abs()/(solution.array().max(1))
).matrix().template lpNorm<p>();
}
//! \brief Project variables on the feasible set
void projection(Eigen::VectorXd& x);
//! \brief Return the step corrected step length if it is too long
double is_step_too_long(Eigen::VectorXd& update);
protected:
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
Log In to Comment