Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F101943472
simplex_program.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, Feb 15, 10:34
Size
4 KB
Mime Type
text/x-c++
Expires
Mon, Feb 17, 10:34 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
24250483
Attached To
rSPECMICP SpecMiCP / ReactMiCP
simplex_program.hpp
View Options
/* =============================================================================
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. *
============================================================================= */
//! \file simplex/simplex_program.hpp
//! \brief Base class for a simplex program
#ifndef SPECMICP_COMMON_SIMPLEX_PROGRAM_HPP
#define SPECMICP_COMMON_SIMPLEX_PROGRAM_HPP
#include "specmicp_common/types.hpp"
namespace
specmicp
{
namespace
simplex
{
//! \brief Base class for a simplex program
//!
//! Uses the CRTP.
//!
//! \f[ \mathrm{min} c^T \cdot x \f]
//! \f[ \mathrm{s.t.} A\cdot x = b
template
<
class
Derived
>
class
SimplexProgram
{
public
:
using
BasisType
=
Eigen
::
Matrix
<
index_t
,
Eigen
::
Dynamic
,
1
>
;
// always the case for now
using
MatrixType
=
Matrix
;
//! \brief Number of equalities constraints
index_t
nb_equalities
()
const
{
return
derived
()
->
nb_equalities_impl
();
}
//! \brief Number of variables
index_t
nb_variables
()
const
{
return
derived
()
->
nb_variables_impl
();
}
//! \brief Number of non zero variables
index_t
size_B
()
const
{
return
nb_equalities
();
}
//! \brief Number of zero variables
index_t
size_N
()
const
{
return
nb_variables
()
-
size_B
();
}
//! \brief Fill the starting guess
//!
//! As B assumed to be the identity value, set the value to b
void
get_starting_guess
(
Eigen
::
Ref
<
Vector
>
m_xB
)
{
derived
()
->
get_starting_guess_impl
(
m_xB
);
}
//! \brief Coefficients of the objective functions
void
get_cB
(
Eigen
::
Ref
<
Vector
>
cB
)
{
derived
()
->
get_cB_impl
(
cB
);
}
//! \brief compute the \f$s_N\f$ vector
//!
//! \f[ c_N - N^T \cdot \lambda \f]
void
compute_sN
(
const
Eigen
::
Ref
<
const
BasisType
>&
Nbasis
,
const
Eigen
::
Ref
<
const
Vector
>&
lambda
,
Eigen
::
Ref
<
Vector
>
sN
)
{
derived
()
->
compute_sN_impl
(
Nbasis
,
lambda
,
sN
);
}
//! \brief Return \f$A(:,q)\f$
void
fill_Aq
(
const
Eigen
::
Ref
<
const
BasisType
>&
Nbasis
,
Eigen
::
Ref
<
Vector
>
Aq
,
index_t
q
)
{
derived
()
->
fill_Aq_impl
(
Nbasis
,
Aq
,
q
);
}
//! \brief Return the coefficient of the objective function corresponding to q
scalar_t
cN_q
(
const
Eigen
::
Ref
<
const
BasisType
>&
Nbasis
,
index_t
q
)
{
return
derived
()
->
cN_q_impl
(
Nbasis
,
q
);
}
private
:
//! \brief Return the derived class
Derived
*
derived
()
{
return
static_cast
<
Derived
*>
(
this
);
}
//! \brief Return a const pointer to the derived class
const
Derived
*
const
derived
()
const
{
return
static_cast
<
const
Derived
*
const
>
(
this
);
}
};
}
}
#endif
// SPECMICP_COMMON_SIMPLEX_PROGRAM_HPP
Event Timeline
Log In to Comment