Algorithms
-
namespace Algorithms
This namespace contains the definitions for the algorithms. For each game, there is a lower-level namespace under which the algorithms are nested.
-
namespace EPEC
Abstact type for other algorithms.
-
class CombinatorialPNE : public Algorithms::EPEC::PolyBase
- #include <epec_combPNE.h>
This class is responsible for the Combinatorial Pure-nash Equilibrium algorithm. In short, it tries to assign to each LCP related to player in the Game::EPEC a single polyhedron. Hence, the resulting Game::NashGame should be easier to solver. If this approximated problem has a solution, it may be a pure-nash equilibrium, since every players’ strategy lays in only one polyhedron.
Public Functions
-
inline virtual void solve()
A general method to solve problems.
-
void solveWithExcluded(const std::vector<std::set<unsigned long int>> &excludeList = {})
Solves the Game::EPEC instance with the algorithm by excluding some combinations of polyhedra that may have been already tested.
- Parameters
excludeList – A set containing the combinations that should be excluded
Private Functions
-
void combPNE(std::vector<long int> combination, const std::vector<std::set<unsigned long int>> &excludeList)
This method initializes the algorithm recursion with
combination
. Each element is the index of a polyhedron for the corresponding player. If the index is -1, then the recursion will generate children for any polyhedron of the given player. Otherwise, if there exist a positive value \(v\) in a location \(l\), player \(l\) will only play strategies in the polyhedron \(v\).- Parameters
combination – A set of either -1 or positive numbers corresponding to the polyhedron of each player. -1 will recurse
excludeList – A set of combinations of polyhedra that should be excluded from the search.
-
inline virtual void solve()
-
class OuterTree
- #include <epec_cutandplay.h>
This class manages the outer approximation tree.
Public Functions
-
inline OuterTree(unsigned int encSize, GRBEnv *env)
Constructor of the Tree given the encoding size.
-
inline void resetFeasibility()
Reset the feasibility parameters for the tree.
-
inline bool getPure() const
Read-only getter for OuterTree::isPure.
-
inline void setFeasible()
Read-only getter for OuterTree::isFeasible.
-
inline void setPure()
Setter for OuterTree::isPure.
-
inline unsigned int getEncodingSize() const
Getter for the encoding size.
-
inline const arma::sp_mat *getV()
Getter for OuterTree::V.
-
inline const arma::sp_mat *getR()
Getter for OuterTree::R.
-
inline unsigned int getVertexCount() const
Getter for OuterTree::VertexCounter.
-
inline unsigned int getRayCount() const
Getter for OuterTree::RayCounter.
-
inline bool addVertex(const arma::vec &vertex, bool checkDuplicates)
Adds a vertex to OuterTree::V.
- Parameters
vertex – The vector containing the vertex
checkDuplicates – True if the method has to check for duplicates
- Returns
True if the vertex was added
-
inline void addRay(const arma::vec &ray)
Adds a ray to OuterTree::R.
- Parameters
ray – The vector containing the ray
-
inline Node *getRoot()
Getter for the root node.
-
inline std::vector<Node> *getNodes()
A getter method for the nodes in the tree.
- Returns
The pointer to the nodes in the tree
-
void denyBranchingLocation(Node &node, const unsigned int &location) const
If a complementarity equation
location
has proven to be infeasible or it isn’t a candidate for branching, this method prevents any further branching on it for the nodenode
.- Parameters
node – The node pointer
location – The denied branching location
-
std::vector<long int> singleBranch(unsigned int idComp, Node &t)
Given the
idComp
and the parent nodet
, creates a single child by branching onidComp
.- Parameters
idComp – The branching id for the complementarity
t – The pointer to the node
- Returns
The node counter stored in a single-element vector
Protected Attributes
-
std::unique_ptr<GRBModel> MembershipLP
Stores the pointer to the MembershipLP associated to the tree.
-
std::unique_ptr<MathOpt::PolyLCP> OriginalLCP
Stores the original LCP. This object is separated from the working one to avoid bugs and numerical problems.
-
arma::sp_mat V = {}
Stores the known extreme vertices of the player’s feasible region. These are used to derive valid cuts, or certify that an equilibrium is inside (outside) the convex-hull of the feasible region.
-
arma::sp_mat R = {}
As in V, but instead of vertices, this object contains rays.
-
unsigned int VertexCounter = 0
The counter for node ids.
-
unsigned int RayCounter = 0
The counter for node ids.
-
bool containsOrigin = false
True if the origin is feasible.
Private Functions
-
inline unsigned int nextIdentifier()
Increments the node counter and get the id of the new node.
Private Members
-
unsigned int EncodingSize = 0
The size of the encoding, namely the number of complementarity equations.
-
unsigned int NodeCounter = 1
The counter for node ids.
-
std::vector<Node> Nodes = {}
Storage of nodes in the tree with the vertices in V
-
bool isPure = {false}
True if the strategy at the current node is a pure-strategy.
-
bool isFeasible{false}
True if the strategy at the current node is feasible for the original game.
-
struct Node
- #include <epec_cutandplay.h>
Public Functions
-
explicit Node(unsigned int encSize)
Constructor for the root node, given the encoding size, namely the number of complementarity equations.
- Parameters
encSize – The number of complementarities
-
Node(Node &parent, unsigned int idComp, unsigned long int id)
Given the parent node address
parent
, theidComp
to branch on, and theid
, creates a new node.- Parameters
parent – The parent node
idComp – The id of the node
id – The The branching candidate
-
Node(Node &parent, std::vector<int> idComps, unsigned long int id)
Given the parent node address
parent
, theidsComp
to branch on (containing all the complementarities ids), and theid
, creates a new node.- Parameters
parent – The parent node pointer
idsComp – The vector of branching locations
id – The node id for the children
-
inline unsigned long int getCumulativeBranches() const
Returns the number of variables that cannot be candidate for the branching decisions, namely the ones on which a branching decision has already been taken, or for which the resulting child node is infeasible.
- Returns
The number of unsuitable branching candidates
-
inline std::vector<bool> getEncoding() const
Getter method for the encoding.
-
inline std::vector<bool> getAllowedBranchings() const
Getter method for the allowed branchings.
Private Members
-
std::vector<unsigned int> IdComps
Contains the branching decisions taken at the node.
-
std::vector<bool> Encoding
An encoding of bool. True if the complementarity condition is included in the current node outer approximation, false otherwise.
-
std::vector<bool> AllowedBranchings
A vector where true means that the corresponding complementarity is a candidate for branching at the current node
-
unsigned long int Id
A long int giving the numerical identifier for the node.
-
Node *Parent = {}
A pointer to the parent node.
Friends
- friend class OuterTree
-
explicit Node(unsigned int encSize)
-
inline OuterTree(unsigned int encSize, GRBEnv *env)
-
class CutAndPlay : public Algorithms::EPEC::PolyBase
- #include <epec_cutandplay.h>
This class is responsible for the Cut-and-Play Algorithm.
Public Functions
-
inline explicit CutAndPlay(GRBEnv *env, Game::EPEC *EPECObject)
Standard constructor.
- Parameters
env – Pointer to the Gurobi environment
EPECObject – Pointer to the EPEC
-
CutAndPlay() = delete
-
inline double getTol() const
Read-Only getter for CutAndPlay::Tolerance.
-
inline void setTol(double tol)
Setter for CutAndPlay::Tolerance.
-
virtual void solve() override
Given the Game::EPEC instance, this method solves it through the outer approximation scheme.
-
void printCurrentApprox()
Prints a log message containing the encoding at the current outer approximation iteration.
-
virtual bool isSolved(double tol = 1e-4) override
Overrides Algorithms::EPEC::PolyBase::isSolved with a custom method.
- Parameters
tol – Numerical tolerance. Currently not useful
- Returns
True if the current outer approximation solution is feasible (then, it is solved)
-
bool isFeasible(bool &addedCuts)
Checks whether the current outer approximation equilibrium is feasible and solves the problem. Otherwise, it adds cuts or generate useful information for the next iterations.
- Parameters
addedCuts – [out] is true if at least a cut has been added
- Returns
-
bool isPureStrategy(double tol = 1e-4) const
Checks whether the current solution is a pure-strategy nash equilibrium.
- Parameters
tol – A numerical tolerance. Currently not used
- Returns
True if the strategy is a pure nash equilibrium
Public Static Functions
-
static void printBranchingLog(std::vector<int> vector)
Given the vector of branching candidates from Algorithms::EPEC::CutAndPlay::getNextBranchLocation, prints a sum up of them.
- Parameters
vector – Output of Algorithms::EPEC::CutAndPlay::getNextBranchLocation
Protected Functions
-
void after()
Standard method for post-solve execution.
-
void updateMembership(const unsigned int &player, const arma::vec &xOfI)
Updates the membership linear-program in the relative Algorithms::EPEC::CutAndPlay::Trees for the player
player
.- Parameters
player – The player index
xOfI – The point for which the membership LP should be updated for
-
int hybridBranching(unsigned int player, OuterTree::Node *node)
Given
player
— containing the id of the player, returns the branching decision for that node given by a hybrid branching rule. In particular, the method return the complementarity id maximizing a combination of constraint violations and number of violated constraints.node
contains the tree’s node. It isn’t const since a branching candidate can be pruned if infeasibility is detected. Note that if the problem is infeasible, namely one complementarity branching candidate results in an infeasible relaxation, then all branching candidates are removed from the list of branching candidates.- Parameters
player – The player id
node – The pointer to the incumbent OuterTree::Node
- Returns
The branching candidate. -1 if none. -2 if infeasible.
-
int infeasibleBranching(unsigned int player, const OuterTree::Node *node)
Given
player
— containing the id of the player, returns the branching decision for that node, where the complementarity is the most (possibly) infeasible one (with both x and z positive). In particular, the method return the (positive) id of the complementarity equation if there is a feasible branching decision atnode
, and a negative value otherwise.- Parameters
player – The player id
node – The pointer to the incumbent OuterTree::Node
- Returns
The branching candidate. Negative if none
-
int deviationBranching(unsigned int player, const OuterTree::Node *node)
Given
player
— containing the id of the player, returns the branching decision for that node, where the complementarity helps include the deviation. In particular, the method return the (positive) id of the complementarity equation if there is a feasible branching decision atnode
, and a negative value otherwise.- Parameters
player – The player id
node – The pointer to the incumbent OuterTree::Node
- Returns
The branching candidate. Negative if none
-
std::unique_ptr<GRBModel> getFeasibilityQP(const unsigned int player, const arma::vec &x)
Given the player index
player
, gets a feasibility quadratic problem enforcingx
to be in the feasible (approximated) region of the Game::EPEC::PlayersQP.- Parameters
player – The player index
x – The strategy for the player
- Returns
A Gurobi pointer to the model
-
void addValueCut(unsigned int player, double RHS, const arma::vec &xMinusI)
Adds a value cut to
player
MathOpt::LCP.- Parameters
player – The index of the player
RHS – The RHS of the value cut
xMinusI – The strategies of the other players
-
bool equilibriumOracle(arma::vec &xOfI, arma::vec &x, unsigned int player, int budget, bool &addedCuts)
The main Equilibrium CutAndPlay loop. Given a player, a maximum number of iterations, a strategy and the other players strategy, it tries to determine if
xOfI
is feasible forplayer
.- Parameters
xOfI – The incumbent strategy for
player
x – The full solution vector
player – The player id
budget – A bound on the number of iteration
addedCuts – The number of added cuts
- Returns
1 if feasible. 0 if infeasible. 2 if iteration limit was hit.
-
bool isFeasiblePure(const unsigned int player, const arma::vec &x)
Given the player index
player
, gets a feasibility quadratic problem enforcingx
to be in the feasible region of the given player.- Parameters
player – The player index
x – The strategy for the player
- Returns
A Gurobi pointer to the model
-
void originFeasibility(unsigned int player)
Gets the LCP for player
player
, and tries to see whether the origin is feasible. If the answer is yes, sets the corresponding object in the tree to true.- Parameters
player – The player’s id
Private Functions
-
std::vector<int> getNextBranchLocation(unsigned int player, OuterTree::Node *node)
Given
player
— containing the id of the player — andnode
containing a node, returns the branching decision for that node, with respect to the current node. In particular, the method return the (positive) id of the complementarity equation if there is a feasible branching decision atnode
, and a negative value otherwise.- Parameters
player – The player id
node – The pointer to the incumbent OuterTree::Node
- Returns
A vector of 4 integers with the branching location given by the most Algorithms::EPEC::CutAndPlay::infeasibleBranching, Algorithms::EPEC::CutAndPlay::deviationBranching, Algorithms::EPEC::CutAndPlay::hybridBranching, and Algorithms::EPEC::CutAndPlay::getFirstBranchLocation, respectively. If an int is negative, there is no real candidate.
-
int getFirstBranchLocation(const unsigned int player, OuterTree::Node *node)
Given
player
— containing the id of the player, returns the branching decision for that node, with no complementarity condition enforced. In particular, the method return the (positive) id of the complementarity equation if there is a feasible branching decision atnode
, and a negative value otherwise.- Parameters
player – The player id
node – The pointer to the incumbent OuterTree::Node
- Returns
The branching candidate. Negative if none
Private Members
-
bool Feasible = {false}
True if a feasible solution has been found.
-
double Tolerance = 3 * 1e-5
The numerical tolerance for the algorithm.
-
inline explicit CutAndPlay(GRBEnv *env, Game::EPEC *EPECObject)
-
class FullEnumeration : public Algorithms::EPEC::PolyBase
- #include <epec_fullenum.h>
This class manages the full enumeration algorithm for Game::EPEC objects.
Public Functions
-
inline FullEnumeration(GRBEnv *env, Game::EPEC *EPECObject)
Standard constructor.
- Parameters
env – Pointer to the Gurobi environment
EPECObject – Pointer to the EPEC
-
virtual void solve()
Solves the Game::EPEC by full enumeration.
-
inline FullEnumeration(GRBEnv *env, Game::EPEC *EPECObject)
-
class InnerApproximation : public Algorithms::EPEC::PolyBase
- #include <epec_innerapp.h>
This class manages the inner enumeration algorithm for Game::EPEC objects. Since each player’s feasible region is a MathOpt::PolyLCP with finitely many polyhedra, each of these region is increasingly expanded with this algorithm. The expansion happens either by adding polyhedra containing profitable moves, or by adding random polyhedra.
Public Functions
-
inline InnerApproximation(GRBEnv *env, Game::EPEC *EPECObject)
Standard constructor.
- Parameters
env – Pointer to the Gurobi environment
EPECObject – Pointer to the EPEC
-
virtual void solve()
A general method to solve problems.
Private Functions
-
void start()
Private main component of the algorithm. Starting from some profitable deviations from an all-zero strategy vector, the algorithm computes the polyhedra containing such deviations, and add them to the approximation. If an approximate equilibrium is found, then the algorithms keeps adding polyhedra by profitable deviation. Otherwise, it adds a random number of Data::EPEC::DataObject::Aggressiveness polyhedra with the method Data::EPEC::DataObject::PolyhedraStrategy.
-
bool addRandomPoly2All(unsigned int aggressiveLevel = 1, bool stopOnSingleInfeasibility = false)
Makes a call to to MathOpt::PolyLCP::addAPoly for each player, and tries to add a polyhedron to get a better inner approximation for the LCP.
aggressiveLevel
is the maximum number of polyhedra it will try to add to each player. Setting it to an arbitrarily high value will mimic complete enumeration.- Parameters
aggressiveLevel – The maximum number of polyhedra to be added to each player
stopOnSingleInfeasibility – If set to true, the function will return false if it cannot add a single polyhedron to a country
- Returns
True when at least a polyhedron is added
-
bool getAllDeviations(std::vector<arma::vec> &deviations, const arma::vec &guessSol, const std::vector<arma::vec> &prevDev = {}) const
Given a potential solution vector
guessSol
, it returns the profitable deviations (if any) for all players indeviations
.- Parameters
deviations – [out] The vector of deviations for all players
guessSol – [in] The guessed solution
prevDev – [in] The previous vector of deviations, if any exist.
- Returns
-
unsigned int addDeviatedPolyhedron(const std::vector<arma::vec> &deviations, bool &infeasCheck) const
Given a vevtor of profitable deviations for all the players, it adds their corresponding polyhedra to the current approximation.
- Parameters
deviations – A vector of vectors containing the deviations
infeasCheck – [out] If at least one player cannot add a polyhedron, the method places false in this output parameter
- Returns
The number of added polyhedra
-
inline InnerApproximation(GRBEnv *env, Game::EPEC *EPECObject)
-
class PolyBase
- #include <epec_polybase.h>
Subclassed by Algorithms::EPEC::CombinatorialPNE, Algorithms::EPEC::CutAndPlay, Algorithms::EPEC::FullEnumeration, Algorithms::EPEC::InnerApproximation
Public Functions
-
inline PolyBase(GRBEnv *env, Game::EPEC *EPECObject)
The standard constructor for a PolyBase algorithm. It creates local MathOpt::PolyLCP objects to work with.
- Parameters
env – The pointer to the Gurobi environment
EPECObject – The pointer to the Game::EPEC object
-
virtual void solve() = 0
A general method to solve problems.
-
bool isSolved(unsigned int *player, arma::vec *profitableDeviation, double tol = -1e-5) const
Checks if Game::EPEC is solved, otherwise it returns a proof.
Analogous to Game::NashGame::isSolved but checks if the given Game::EPEC is solved. If it is, then returns true. If not, it returns the country which has a profitable deviation in
player
and the profitable deviation inprofitableDeviation
.Tolerance
is the tolerance for the check. If the improved objective after the deviation is less thanTolerance
, then it is not considered as a profitable deviation.Thus we check if the given point is an \(\epsilon\)-equilibrium. Value of \(\epsilon \) can be chosen sufficiently close to 0.
Warning
Setting
Tolerance
= 0 might even reject a real solution as not solved. This is due to Numerical issues arising from the LCP solver (Gurobi).- Parameters
player – The id of the player
profitableDeviation – An output (possibly non-changed) vector containing the profitable deviation for the given player
tol – A numerical tolerance
- Returns
True if there is no profitable deviation, namely the player is optimal
-
virtual bool isSolved(double tol = 1e-5)
A method to check whether the EPEC is solved or not, given a numerical tolerance.
Checks whether the current Game::EPEC instance is solved for any player, up to a numerical tolerance.
- Parameters
tol – The numerical tolerance
- Returns
True if the game is solved
-
void makeThePureLCP()
Creates an LCP for the inner-full approximation schemes of MathOpt::PolyLCP that explicitly searches for pure-strategy equilibria. The original LCP is moved to Game::EPEC::LCPModelBase.
-
double getValLeadFollPoly(unsigned int i, unsigned int j, unsigned int k, double tol = 1e-5) const
For the
i
-th leader, gets thek
-th pure strategy at positionj
.- Parameters
i – The leader index
j – The position index
k – The pure strategy index
tol – A numerical tolerance
- Returns
The queried attribute
-
double getValLeadLeadPoly(unsigned int i, unsigned int j, unsigned int k, double tol = 1e-5) const
For the
i
-th leader, gets thek
-th pure strategy at leader positionj
.- Parameters
i – The leader index
j – The position index
k – The pure strategy index
tol – A numerical tolerance
- Returns
The queried attribute
-
double getValProbab(unsigned int i, unsigned int k) const
The probability associated with the
k
-th polyhedron of thei
-th leader.- Parameters
i – The index of the player
k – The index of the polyhedron
- Returns
The queried attribute
-
bool isPureStrategy(unsigned int i, double tol = 1e-5) const
Checks whether the current strategy for the
i
player is a pure strategy.- Parameters
i – The index of the player
tol – A numerical tolerance
- Returns
True if it is a pure equilibrium strategy
-
bool isPureStrategy(double tol = 1e-5) const
Checks if the current equilibrium strategy in Game::EPEC is a pure strategy.
- Parameters
tol – A numerical tolerance
- Returns
True if it is a pure equilibrium strategy
-
std::vector<unsigned int> mixedStrategyPoly(unsigned int i, double tol = 1e-5) const
Returns the indices of polyhedra feasible for the leader, from which strategies are played with probability greater than the tolerance.
- Parameters
i – The index of the player
tol – A numerical tolerance
- Returns
The indices of polyhedra with active probabilities
-
unsigned int getPositionLeadFollPoly(unsigned int i, unsigned int j, unsigned int k) const
Get the position of the
k
-th follower variable of thei
-th leader, in thej
-th feasible polyhedron.- Parameters
i – The leader index
j – The polyhedron index
k – The follower variable index
- Returns
The position for the queried attribute
-
unsigned int getPositionLeadLeadPoly(unsigned int i, unsigned int j, unsigned int k) const
Get the position of the
k
-th leader variable of thei
leader, in thej-th
feasible polyhedron.- Parameters
i – The leader index
j – The polyhedron index
k – The leader variable index
- Returns
The position for the queried attribute
-
unsigned long int getNumPolyLead(unsigned int i) const
Get the number of polyhedra used in the approximation for the
i
leader.- Parameters
i – The leader index
- Returns
The queried number of polyhedra
-
unsigned int getPositionProbab(unsigned int i, unsigned int k) const
Get the position of the probability associated with the
k
-th polyhedron (k
-th pure strategy) of thei
-th leader. However, if the leader has an inner approximation with exactly 1 polyhedron, it returns 0;.- Parameters
i – The leader index
k – The polyhedron index
- Returns
The probability position associated with the queried values
Protected Functions
-
inline void after()
This method is called after the PolyBase::solve operation. It fills statistics and can be forcefully overridden by inheritors. The responsibility for calling this method is left to the inheritor.
Protected Attributes
-
GRBEnv *Env
This is the abstract class manages the algorithms for Game::EPEC. Since they are all based on MathOpt::PolyLCP, the class keeps a local copy of objects of that class. It provides a constructor where the Gurobi environment and the EPEC are passed.
A pointer to the Gurobi Environment
-
Game::EPEC *EPECObject
A pointer to the Game::EPEC instance.
-
std::vector<std::shared_ptr<MathOpt::PolyLCP>> PolyLCP = {}
Local MathOpt::PolyLCP objects.
-
inline PolyBase(GRBEnv *env, Game::EPEC *EPECObject)
-
class CombinatorialPNE : public Algorithms::EPEC::PolyBase
-
namespace IPG
-
class Algorithm
- #include <ipg_algorithms.h>
Subclassed by Algorithms::IPG::CutAndPlay, Algorithms::IPG::ZERORegrets
Public Functions
-
virtual void solve() = 0
A method to solve IPGs.
-
virtual bool isSolved() const = 0
A method to check whether the IPG is solved or not, given a numerical tolerance.
-
virtual bool isPureStrategy() const = 0
A method to check whether the IPG solution is a pure equilibrium or not, given a numerical tolerance.
-
inline double getTol() const
-
inline void setTol(double tol)
Protected Attributes
-
GRBEnv *Env
-
bool Solved = {false}
True if the IPG has been solved.
-
bool Pure = {false}
True if all the players are playing a pure strategy.
-
bool Infeasible = {false}
True if the game is infeasible.
-
double Tolerance = 1e-6
The numeric tolerance.
-
virtual void solve() = 0
-
struct IPG_Player
- #include <ipg_cutandplay.h>
Public Functions
-
~IPG_Player() = default
-
inline IPG_Player(unsigned int incumbentSize, double &tol)
-
bool addVertex(const arma::vec &vertex, const bool checkDuplicate = false)
Given
vertex
, it adds a vertex to the field V. IfcheckDuplicate
is true, it will check whether the vertex is already contained in the pool.- Parameters
vertex – The input vertex
checkDuplicate – True if the method needs to check for duplicates
- Returns
True if the vertex is added.
-
bool addRay(const arma::vec &ray, const bool checkDuplicate = false)
Given
ray
, it adds a ray to the field R. IfcheckDuplicate
is true, it will check whether the ray is already contained in the pool.- Parameters
ray – The input ray
checkDuplicate – True if the method needs to check for duplicates
- Returns
True if the ray is added.
-
bool addCuts(const arma::sp_mat &LHS, const arma::vec &RHS)
Given
LHS
,RHS
, it adds the inequalities to the field CutPool_A and b, the CoinModel, and the working IP_Param.- Parameters
LHS – The LHS matrix
RHS – The RHS vector
- Returns
true if the inequality is added.
Protected Attributes
-
std::unique_ptr<GRBModel> MembershipLP = {}
This structure manages the IPG data for each player of the game, given the CutAndPlay.
The model approximating the feasible region with vertices and rays
-
std::shared_ptr<MathOpt::IP_Param> ParametrizedIP = {}
The (working) player integer program, to which cuts are added.
-
std::shared_ptr<OsiGrbSolverInterface> CoinModel = {}
Quick workaround for now. This object stores the CoinOR model related to the field ParametrizedIP.
-
arma::sp_mat V = {}
This object stores an array of points — for each player — that are descriptor for the convex-hull of the integer programming game.
-
arma::sp_mat R = {}
As in V, but for rays.
-
bool containsOrigin = false
True if the origin is a feasible point.
-
unsigned int VertexCounter = 0
The number of Vertices in the membership LP.
-
unsigned int RayCounter = 0
The number or Rays in the membership LP.
-
arma::sp_mat CutPool_A = {}
Stores the LHS of the valids cuts for the convex hull of the player’s IPG.
-
arma::vec CutPool_b = {}
Stores the RHS of the valids cuts for the convex hull of the player’s IPG.
-
double Tolerance = 1e-6
Numerical tolerance.
-
arma::vec Incumbent
Stores the current strategy of the player at a given iteration.
-
arma::vec DualIncumbent
Stores the (dual) current strategy of the player at a given iteration.
-
double Payoff
Stores the current payof.
-
bool Pure
True if the strategy is pure.
-
bool Feasible = false
Friends
- friend class Algorithms::IPG::CutAndPlay
-
~IPG_Player() = default
-
class CutAndPlay : public Algorithms::IPG::Algorithm
- #include <ipg_cutandplay.h>
This class is responsible for the Cut-and-Play algorithm for IPG.
Public Functions
-
inline CutAndPlay(GRBEnv *env, Game::IPG *IPGObj)
Standard constructor.
- Parameters
env – The Gurobi environment
IPGObj – The IPG object
-
virtual void solve()
Solves the IPG with the Equilibrium CutAndPlay algorithm.
-
inline virtual bool isSolved() const
A method to check whether the IPG is solved or not, given a numerical tolerance.
-
virtual bool isPureStrategy() const
Returns true if all players are playing a pure strategy in a Nash Equilibrium.
- Returns
True if the Equilibrium is pure
Private Functions
-
void initialize()
The last Nash Game to which the LCP object is associated.
This method initializes some fields for the algorithm. Also, it warm starts the initial strategies to pure best responses.
-
arma::vec buildXminusI(const unsigned int i)
Given the player id
i
, builds the vector x^{-i} from the current working strategies.- Parameters
i – The player id
- Returns
The other players strategies (except
i
)
-
void initializeEducatedGuesses()
Initializes some pure-strategies for each player.
-
void initializeCoinModel(const unsigned int player)
This method builds the Coin-OR model used in CutAndPlay::externalCutGenerator for the given player.
- Parameters
player – The player’s id
-
unsigned int externalCutGenerator(unsigned int player, int maxCuts, bool rootNode, bool cutOff)
Given a player
player
, a number of maximum cuts to generatemaxcuts
and a boolrootNode
, this method generates some valid inequalities for theplayer
‘s integer program. This method uses Coin-OR CGL. So far, Knapsack covers, GMI and MIR inequalities are used. Also, note that there is no recursive cut generation (meaning, we do not generate cuts from a previous cut pool) as to better manage numerical stability.cutOff
requires to cut off the current solution forplayer
.- Parameters
player – The current player id
maxCuts – The maximum number of cuts
rootNode – True if the cut generation happens at the root node
cutOff – True if the cuts are required to cutoff the current solution
- Returns
The number of added cuts
-
bool addValueCut(unsigned int player, double RHS, const arma::vec &xMinusI)
Given a player
player
, one of its best responsesxOfIBestResponses
, the strategies of the other playersxMinusI
, it adds an inequality of the type.\[ f^i(x^i,\bar x^{-i}) \geq f^i(\hat x^i, \bar x^{-i})\]to the cut pool of that player.- Parameters
player – The player id
RHS – The RHS value
xMinusI – The input
\[ x^{-i} \]
- Returns
True if the cut was added
-
int preEquilibriumOracle(const unsigned int player, int &addedCuts, arma::vec &xOfI, arma::vec &xMinusI)
Given the player id
player
, checks whether the current strategy is feasible or not. In order to do so, a more complex separation technique may be called.- Parameters
player – The player id
addedCuts – Filled with the number of added cuts
xOfI – The strategy of
player
xMinusI – The strategy of the other players
- Returns
1 if feasible. 0 if infeasible. 2 if iteration limit was hit.
-
void updateMembership(const unsigned int &player, const arma::vec &vertex)
Update the Membership Linear Program for the given player and the verter
vertex
.- Parameters
player – The player id
vertex – The input point to be checked
-
int equilibriumOracle(const unsigned int player, const unsigned int iterations, const arma::vec &xOfI, const arma::vec &xMinusI, int &addedCuts)
The main Equilibrium CutAndPlay loop. Given a player, a maximum number of iterations, a strategy and the other players strategy, it tries to determine if
xOfI
is feasible forplayer
.- Parameters
player – The player id
iterations – The bound on iterations
xOfI – The strategy of
player
xMinusI – The strategies of other players
addedCuts – The number of added cuts
- Returns
1 if feasible. 0 if infeasible. 2 if iteration limit was hit.
-
bool checkTime(double &remaining) const
Checks if there is more time remaining.
- Parameters
remaining – An output filled with the time remaining
- Returns
True if there is still time left.
-
void initLCPObjective()
Initialize the LCP Objective with the quadratic and linear terms. These will be later used if necessary.
-
ZEROStatus equilibriumLCP(double localTimeLimit, bool build = true, bool firstSolution = true)
Creates and solves the equilibrium LCP wrt the current game approximation.
- Parameters
localTimeLimit – A time limit for the computation
build – If true, a new LCP will be built. Otherwise, the last one will be used.
firstSolution – If true, the method will just seek for one solution.
- Returns
The ZEROStatus for the equilibrium LCP
Private Members
-
arma::sp_mat LCP_Q
Quadratic matrix for the LCP objective.
-
arma::vec LCP_c
Linear vector for the LCP objective.
-
std::vector<std::unique_ptr<IPG_Player>> Players
The support structures of IPG_Players.
-
std::vector<std::pair<std::string, int>> Cuts
Log of used cutting planes.
-
arma::vec zLast
The last z solution. Useful for warmstarts.
-
arma::vec xLast
The last x solution. Useful for warmstarts.
-
double objLast = -GRB_INFINITY
Last objective from the equilibrium LCP. Used as cutOff.
Friends
- friend class Game::IPG
-
inline CutAndPlay(GRBEnv *env, Game::IPG *IPGObj)
-
class ZERORegrets : public Algorithms::IPG::Algorithm
- #include <ipg_zeroregrets.h>
This class is responsible for the ZERORegrets algorithm (Dragotto and Scatamacchia) for IPGs.
Public Functions
-
inline ZERORegrets(GRBEnv *env, Game::IPG *IPGObj)
Standard constructor.
- Parameters
env – The Gurobi environment
IPGObj – The IPG object
-
virtual void solve()
Solves the IPG with the Equilibrium CutAndPlay algorithm.
-
inline virtual bool isSolved() const
A method to check whether the IPG is solved or not, given a numerical tolerance.
-
inline virtual bool isPureStrategy() const override
A method to check whether the IPG solution is a pure equilibrium or not, given a numerical tolerance.
Private Functions
-
void initialize()
This method initializes some fields for the algorithm. Also, it warm starts the initial strategies to pure best responses.
-
bool addEquilibriumInequality(unsigned int player, const arma::vec &xOfI)
Generates an equilibrium inequality starting from xOfI.
- Parameters
player –
xOfI –
- Returns
-
bool checkTime(double &remaining) const
Checks if there is more time remaining.
- Parameters
remaining – An output filled with the time remaining
- Returns
True if there is still time left.
-
ZEROStatus equilibriumMIP(double localTimeLimit)
Private Members
-
std::unique_ptr<GRBModel> JointProgram = {}
The joint MIP program to which cuts are added.
-
std::vector<std::pair<std::string, int>> Cuts
Log of used cutting planes.
-
std::vector<arma::vec> xLast
The last x solution, for each player.
-
GRBVar **x = {}
-
GRBVar *p = {}
-
double objLast = +GRB_INFINITY
Last objective from the equilibrium MIP. Used as cutOff.
Friends
- friend class Game::IPG
-
inline ZERORegrets(GRBEnv *env, Game::IPG *IPGObj)
-
class Algorithm
-
namespace EPEC