Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F91279815
micpsolver.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
Sat, Nov 9, 15:10
Size
10 KB
Mime Type
text/x-c++
Expires
Mon, Nov 11, 15:10 (2 d)
Engine
blob
Format
Raw Data
Handle
22231457
Attached To
rSPECMICP SpecMiCP / ReactMiCP
micpsolver.hpp
View Options
/* =============================================================================
Copyright (c) 2014 - 2016
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 SPECMIC_MICPSOLVER_MICPSOLVER_HPP
#define SPECMIC_MICPSOLVER_MICPSOLVER_HPP
#include "specmicp_common/types.hpp"
#include "micpsolver_base.hpp"
//! \file micpsolver.hpp
//! \brief The MiCP solver
namespace
specmicp
{
//! \namespace specmicp::micpsolver
//! \brief the MiCP Solver(s) and related classes
namespace
micpsolver
{
//! \enum ReformulationF
//! \brief The different reformulation functions
//!
enum
class
ReformulationF
{
CCK
,
//!< The CCK reformulation for MiCP problems
BoxVI
//!< Special reformulation for box-VI problems
};
//! \class MiCPSolver
//! \brief The MiCP Solver
//!
//! ### Mathematics
//!
//! Solve
//! - \f$ G(u, v) = 0\f$
//! - \f$ 0 \leq v \perp H(u,v) \geq 0 \f$
//!
//! Using the penalized Fisher-Burmeister C-function.
//! This is a semismooth Newton solver with many features including :
//!
//! - Non-monotone linesearch
//! - Scaling of the jacobian
//! - Check of the condition number of the jacobian
//! - Check of the descent directions
//! - Crashing
//! - ...
//!
//! References :
//! - Munson et al (2001) \cite Munson2001
//! - Facchinei and Pang (2003) \cite Facchinei2003
//!
//! ### Usage
//!
//! The options to customize the solver are contained in the struct
//! specmicp::micpsolver::MiCPSolverOptions .
//!
//! \tparam program_t a subclass of MiCPProg
//!
//! To separate the program and the solver, the program is passed as a shared_ptr.
//! The MiCPSolver is implemented using the CRTP pattern.
//! The advantage is that we have compile time polymorphism instead of
//! the runtime polymorphism when using the virtual functions.
//! I am not sure of the effectiveness of this method (I have not run any benchmarks)
//! but since the solver is destined to be called many times any small gain is good to have.
//! Also, more compile-time optimization is possible (especially inlining).
//! Anyway, that was fun to code...
//!
//! \sa specmicp::micpsolver::MiCPProgram
template
<
class
program_t
,
ReformulationF
reform_f
=
ReformulationF
::
CCK
>
class
MiCPSolver
:
public
MiCPSolverBaseProgram
<
program_t
>
{
public
:
using
base
=
MiCPSolverBaseProgram
<
program_t
>
;
using
base
::
get_options
;
using
base
::
get_program
;
using
base
::
get_perfs
;
using
base
::
get_neq
;
using
base
::
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
);
}
// 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
);
//! \brief Setup the jacobian
//!
//! \param[in] x the current solution
void
setup_jacobian
(
Eigen
::
VectorXd
&
x
);
//! \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 specmicp::MiCPSolver::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,out] update the update of the timestep
//! \param[in,out] x the current solution
//!
//! References :
//! - \cite Dennis1983
//! - \cite Munson2001
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
);
//! \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
)
{
Reformulation
<
reform_f
>::
reformulate_residuals
(
this
,
x
,
r
,
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
)
{
Reformulation
<
reform_f
>::
reformulate_residuals_inplace
(
this
,
x
,
r
);
}
//! \brief Reformulation of the jacobian
//!
//! r is the original vector of residuals (not reformulated)
void
reformulate_jacobian
(
const
Vector
&
x
)
{
Reformulation
<
reform_f
>::
reformulate_jacobian
(
this
,
x
);
}
protected
:
//! \brief Structure in charge of the reformulation
template
<
ReformulationF
ref_func
,
int
dummy
=
0
>
struct
Reformulation
{
//! \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
static
void
reformulate_residuals
(
MiCPSolver
*
const
parent
,
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
static
void
reformulate_residuals_inplace
(
MiCPSolver
*
const
parent
,
const
Vector
&
x
,
Vector
&
r
);
//! \brief Reformulation of the jacobian
//!
//! r is the original vector of residuals (not reformulated)
static
void
reformulate_jacobian
(
MiCPSolver
*
const
parent
,
const
Vector
&
x
);
};
//! \brief Return the residuals (as givent by the system)
Vector
&
get_residuals
()
{
return
m_residuals
;
}
//! \brief Return the reformulated residuals
Vector
&
get_reformulated_residuals
()
{
return
m_phi_residuals
;
}
//! \brief Return the jacobian
//!
//! The reformulation is in-place so the meaning is dependent on the call
Matrix
&
get_jacobian
()
{
return
m_jacobian
;
}
private
:
// Residuals and jacobian
Vector
m_residuals
;
//!< The residuals
Vector
m_phi_residuals
;
//!< The reformulated residuals
Vector
m_grad_phi
;
//!< The gradient of the reformulated residuals
Matrix
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
{
false
};
//!< True if the update was computed using the gradient
};
}
// end namespace micpsolver
}
// end namespace specmicp
// ###############//
// Implementation //
// ###############//
#include "reformulations.inl"
// definitions of the reformulations
#include "micpsolver.inl"
// implementations of the solver
#endif
// SPECMIC_MICPSOLVER_MICPSOLVER_HPP
Event Timeline
Log In to Comment