Content-type: text/html; charset=UTF-8 Man page of lrs

lrs

Section: lrs 7.3 (1)
Updated: 2023.1.4
Index Return to Main Contents
 

Name

lrs - Convert between representations of convex polyhedra, remove redundant inequalities, convex hull computation, volume, triangulation, solution to linear programs in exact precision.  

Synopsis

lrs|plrs [input-file] [output-file]
redund|minrep [input-file] [output-file] See: redund(1)
fel [input-file] [output-file] See: fel(1)
hvref/xref [input-file] See: xref(1)
 

Description

These programs are part of and must be compiled with lrslib which is a C library. All computations are done in exact arithmetic.

A polyhedron can be described by a list of inequalities (H-representation) or by a list of its vertices and extreme rays (V-representation).

lrs converts an H-representation of a polyhedron to its V-representation and vice versa, known respectively as the vertex enumeration and facet enumeration problems (see Example (1) below). For V-representations the volume can be computed and a triangulation produced. lrs can also be used to solve a linear program, remove linearities from a system, and extract a subset of columns.
plrs is a primitive multi-threaded version of lrs that gives some speedup on shared memory machines. See threads option. For a fully parallel version see mplrs(1)

redund|minrep
These are aliases of lrs which can also perform these functions via the redund, redund_list and testlin options. For a more complete description and examples see redund(1).
H-representation: redundant inequalities in an input H-representation are removed and the remaining inequalities output. In addition minrep checks for hidden linearities.
V-representation: outputs all extreme points and extreme rays, often called the convex hull problem. In addition minrep checks for hidden linearities.

fel projects an input H-representation onto a given set of variables using Fourier-Motzkin elimination. For a V-representation it extracts the specified columns. The output is non-redundant and can be can be piped directly into lrs. This is an alias of lrs which can also perform these functions via the eliminate and project options.

hvref/xref produce a cross reference list between H- and V-representations. See xref(1).

mplrs is Skip Jordan's parallel wrapper based on MPI for lrs/redund using the same input and output formats. See: mplrs(1)

Fukuda's FAQ page[1] contains a more detailed introduction to the problem, along with many useful tips for the new user. User's guide for lsrslib[8]

 

File formats

File formats were developed jointly with Komei Fukuda and are compatible with cdd/cddlib[2].
The input for lrs/redund is an H- or V-representation of a polyhedron.


 name
 H-representation [or V-representation]
 {options}
 {linearities}
 begin
  m n rational [or integer]
 {input matrix}
 end 
 {options}

name is a user supplied name for the polyhedron.  Comments may appear before the begin or after the end and should begin with a special character such as "*".

If the representation is not specified H-representation is assumed. The input coefficients are read in free format, and are not checked for type. Coefficients are separated by white space. m is the number of rows and n the number of columns of the input matrix.  

H-representation

m is the number of input rows, each being an inequality or equation.
n is the number of input columns and d=n-1 is the dimension of the input.
An inequality or equation of the form:

b + a_1 x_1 + ... + a_d x_d >=  0

b + a_1 x_1 + ... + a_d x_d =  0

is input as the line:

b  a_1 ... a_d

The coefficients can be entered as integers or rationals in the format x/y. To distinguish an equation a linearity option must be supplied before the begin line (see below).  

V-representation

m is the number of input rows, each being a vertex, ray or line.
n is the number of input columns and d=n-1 is dimension of the input.
Each vertex is given in the form:

1   v_1   v_1 ... v_d

Each ray is given in the form:

0   r_1   r_2...   r_d

where r_1  ...   r_d is a point on the ray.

There must be at least one vertex in each file. For bounded polyhedra there will be no rays entered. The coefficients can be entered as integers or rationals in the format x/y. An input line can be specified as a ray and then included in the linearity option (see below).

Note for cdd users: Note the input files for lrs are read in free format. lrs will look for exactly m*n rationals or integers separated by white space (blank, carriage return, tab etc.). lrs will not "drop" extra columns of input if n is less than the number of columns supplied.

 

Options

Almost all options are placed after the end statement, maintaining compatibility with cdd. Where this is not the case, it will be mentioned explicitly.

allbases This option instructs lrs to list each vertex (or facet) for each of its bases. This option is often combined with printcobasis.

bound  x (H-representation only). Either the maximize or minimize option should be selected. x is an integer or rational. For maximization (resp. minimization) the reverse search tree is truncated  whenever the current objective value is less (resp. more) than x.

cache n     (default n=50)
lrs stores the latest  n dictionaries in the reverse search tree. This speeds up the backtracking step, but requires more memory.

debug  startingcobasis endingcobasis
Print out cryptic but detailed trace, dictionaries etc. starting at #B=startingcobasis and ending at #B=endingcobasis. debug 0 0 gives a complete trace.

digits n (lrsmp arithmetic only - placed before the begin statement)
n is the maximum number of decimal digits to be used. If this is exceeded the program terminates with a message and can usually be restarted with the restart option. The default is set to 100 digits. At the end of a run a message is given informing the user of the maximum integer size encountered.

dualperturb If lrs is executed with the maximize or minimize option, the reverse search tree is rooted at an optimum vertex for this function. If there are multiple optimum vertices, the output will often not be complete. This option gives a small perturbation to the objective to avoid this. A warning message is given if the starting dictionary is dual degenerate.

estimates k
Estimate the output size. Used in conjunction with maxdepth and seed. See: Estimation[3].

eliminate k i_1 i_2 ... i_k (see fel(1))
(H-representation) Eliminates k variables in an H-representation corresponding to cols i_1 .. i_k by projection onto the remaining variables using the Fourier-Motzkin method. Variables are eliminated in the order given and redundancy is removed after each iteration.
(V-representation) Delete the k given columns from the input matrix and remove redundancies (cf. extract where redundancies are not removed).
Column indices are between 1 and n-1 and column zero cannot be eliminated. The output is a valid lrs input file. See also project and extract

extract [ k i_1 i_2 ... i_k ]
(H-representation) A preprocessing step to remove linearities (if any) in an H-representation and resize the A matrix. The output as a valid lrs input file. The resulting file will not contain any equations but may not be full dimensional as there may be additional linearities in the remaining inequalities. Options in the input file are stripped. The user can specify the k columns i_1 i_2 ... i_k to retain otherwise if k=0 or is ommitted the lex-min set of columns are retained, equivalent to the the order 1,2,..n-1. Linearly dependent columns are skipped and additional indices are taken from 1,2,...,n-1 as necessary. If there are no linearities in the input file the given columns are retained and the other ones are deleted.
(V-representation) Extract the given columns from the input file outputing a valid lrs input file. Options are stripped.

geometric   (H-representation  or voronoi option only) Each ray is printed together with the vertex with which it is incident.

incidence This option automatically switches on printcobasis. For input H-representation, indices of all input inequalities that contain the vertex/ray that is about to be output. For input V-representation, indices of all input vertices/rays that lie on the facet that is about to be output. A starred index indicates that this vertex  is also in the cobasis, but is not contained in the facet. It arises due to the lifting operation used with input V-representations.

linearity  k  i_1  i_2  ...  i_k (placed before begin statement)
(H-representation) The k rows i_1  i_2  ...  i_k  of the input file represent equations. (V-representation) The k rows, which should have a zero in column 1, represent lines in space (rather than rays).

lponly Solve the LP given by the input H-representation with objective function specified by the maximize or minimize options and terminate. Use with verbose option to get dual variables. See: Linear Programming[4].

maxcobases k
The search will be truncated after k cobases have been generated.

maxdepth k
The search will be truncated at depth k. All bases with depth less than or equal to k will be computed.  k is  a non-negative integer, and this option is used for estimates - see Estimation[3]. Note: For H-representations, rays at depth k will not be reported. For V-representations, facets at depth k will not be reported.

maximize  b  a_1 ... a_{n-1}  
minimize  b  a_1 ... a_{n-1}  
(H-representation) The starting vertex maximizes (or minimizes) the function  b + a_1 x_1+ ... + a_{n-1} x_{n-1}.
The dualperturb option may be needed to avoid dual degeneracy. Often used with lponly.
(V-representation, v.7.2) The input file row numbers maximizing(minimizing) the function are output along with the optimum value. Using verbose the optimizing lines are also printed. With minimization a facet gives an optimum value of zero, a negative value indicates infeasibility and a positive value indicates strong redundancy.

maxincidence n k
Prunes the search tree when the depth is at least k and the current vertex/facet has incidence at least n. Using verbose a message is printed whenever the search tree is pruned.

maxoutput n
Limits number of output lines produced (either vertices+rays or facets) to n

mindepth k
Backtracking will be terminated at depth k.

nonnegative (This option must come before the begin statement - H-representation only)   Bug: Can only be used if the origin is a vertex of the polyhedron  For problems where the input is an H-representation of the form b+Ax>=0, x>=0 (ie. all variables non-negative, all constraints inequalities) it is not necessary to give the non-negative constraints explicitly if the nonnegative option is used. This option cannot be used for V-representations, or with the linearity option (in which case the linearities will be treated as inequalities). This option may be used with redund , but the implied nonnegativity constraints are not tested themselves for redundancy.

project k i_1 i_2 ... i_k (see fel(1))
(H-representation) Project the polyhedron onto the k variables corresponding to cols i_1 .. i_k using the Fourier-Motzkin method. Column indices are between 1 and n-1 and column zero is automatically retained. Variables not contained in the list are eliminated using a heuristic which chooses the column which minimizes the product of the number of positive and negative entries. Redundancy is removed after each iteration using linear programming.
(V-representation) Extract the k given columns from the input matrix and remove redundancies. Column indices are between 1 and n-1 and column zero is automatically extracted (cf. extract where redundancies are not removed).
The output as a valid lrs input file. See also eliminate and extract

printcobasis k
Every k-th cobasis is printed. If k is omitted, the cobasis is printed for each vertex/ray/facet that is output. For a long run it is useful to print the cobasis occasionally so that the program can be restarted if necessary. H-representation: the cobasis is a list the indices of the inequalities from the input file that define the current vertex or ray. For rays the cobasis is the cobasis of the vertex from which the ray emanates. One of the indices is starred, this indicates the inequality to be dropped from the cobasis to define the ray. If the allbases option is used, all cobases will be printed. V-representation: the cobasis is a list of the input vertices/rays that define the current facet. See option incidence for more information.

printslack (H-representation only) A list of the indices of the input inequalities that are satisfied strictly for the current vertex, ie. corresponding slack variable is positive. If nonnegative is set, the list will also include indices n+i for each decision variable x_i which is positive.

redund start end (see redund(1))
Check input lines with line numbers from start to end and remove any redundant lines.
redund 0 0 will check all input lines.

redund_list k i_1 i_2 ... i_k (see redund(1))
Check the k input line numbers with indices i_1 i_2 ... i_k and remove any redundant lines.

restart  V# R# B# depth {facet #s or vertex/ray #s}
lrs can be restarted from any known cobasis. The calculation will proceed to normal termination. All of the information is contained in the output from a printcobasis option.  The order of the indices is very important, enter them exactly as they appear in the output from the previously terminated run.

seed k
Set the random number generator seed=k. Used with estimates.

startingcobasis i_1  i_2  ...  i_{n-1}
lrs will start from the given cobasis which which is a list of the inequalities (for H-representation) or vertices/rays (for V-representation) that define it. If it is invalid, or this option is not specified, lrs will find its own starting cobasis.

testlin (before the begin line only) H-representation only (new 7.3)
redund An LP test will be made for hidden linearities at the beginning of the run and is reported. If there are no hidden linearities one LP per constraint tests for redundancy. (mplrs can be run with no verification step.) If hidden linearities exist two LPs per constraint search for hidden linearities and remove redundancies. In both cases the run ends with a minimum set of linearities and inequalities (ie. no hidden inequalities or duplicates) and the dimension is reported. (mplrs will find all linearities but some non redundant inequalities may be eliminated, so a second run is required.)
lrs If neither redund or redund_list options are present the initial LP test is made, reported and the run halted. Otherwise same as redund above.

threads n (new in 7.3, currently for plrs1,plrs2,plrsgmp,plrsmp and script plrs)
Multithreaded lrs using openMP for a parallel for loop at depth 0 in the search tree. Disabled for mplrs. If n is not specified the default openMP max threads is used.

truncate The reverse search tree is truncated(pruned)  whenever a new vertex is encountered. Note: This does note necessarily produce the set of all vertices adjacent to the optimum vertex in the polyhedron, but just a subset of them.

verbose Print slightly more detailed information about the run.

volume (V-representation only) Compute the volume and, if the verbose option is also included, output a triangulation. See Volume Computation[5].

voronoi (V-representation  only - place immediately after end statement)
Compute Voronoi diagram - see Voronoi Diagrams[6].  

Arithmetic

From version 7.1 lrs/redund/mplrs use hybrid arithmetic with overflow checking, starting in 64bit integers, moving to 128bit (if available) and then GMP. Overflow checking is conservative to improve performance: eg. with 64 bit arithmetic, a*b triggers overflow if either a or b is at least 2^31, and a+b triggers an overflow if either a or b is at least 2^62. Typically problems that can be solved in 64bits run 3-4 times faster than with GMP and inputs solvable in 128bits run twice as fast as GMP.

Various arithmetic versions are available and can be built from the makefile:

lrs1 Fixed length 64 bit integer arithmetic, terminates on overflow.

lrs2 Fixed length 128 bit integer arithmetic, terminates on overflow.

lrsmp Built in extended precision integer arithmetic, uses digits option above.

lrsgmp GNU MP which must be installed first from https://gmplib.org/.

lrsflint FLINT hybrid arithmetic which must be installed first from http://www.flintlib.org/

 

Examples

(1) Convert the H-representation of a cube given cube by 6 the six inequalities
-1 <= x_i <= 1 , i=1,2,3 into its V-representation consisting of 8 vertices.


 % cat cube.ine
 cube.ine
 H-representation
 begin
 6 4 rational
 1  1  0  0
 1  0  1  0
 1  0  0  1
 1 -1  0  0
 1  0  0 -1
 1  0 -1  0
 end


 % lrs cube.ine


 *lrs:lrslib v.6.3 2018.4.11(64bit,lrslong.h,overflow checking)
 *Input taken from file cube.ine
 cube.ine
 V-representation
 begin
 ***** 4 rational
 1  1  1  1
 1 -1  1  1
 1  1 -1  1
 1 -1 -1  1
 1  1  1 -1
 1 -1  1 -1
 1  1 -1 -1
 1 -1 -1 -1
 end
 *Totals: vertices=8 rays=0 bases=8 integer_vertices=8

 

Notes

1.
FAQ page
https://inf.ethz.ch/personal/fukudak/polyfaq/polyfaq.html
2.
cdd
https://inf.ethz.ch/personal/fukudak/cdd_home/
3.
Estimation.
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html#Estimation
4.
Linear Programming
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html#Linear%20Programming
5.
Volume Computation.
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html#Volume%20Computation
6.
Voronoi Diagrams.
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html#Voronoi%20Diagrams
7.
redund: extreme point enumeration and eliminating redundant inequalities
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html#redund
8.
User's guide for lrslib
http://cgm.cs.mcgill.ca/%7Eavis/C/lrslib/USERGUIDE.html
 

Author

David Avis <avis at cs dot mcgill dot ca >  

See also

mplrs(1), lrslib(5), lrsnash(1)


 

Index

Name
Synopsis
Description
File formats
H-representation
V-representation
Options
Arithmetic
Examples
Notes
Author
See also

This document was created by man2html, using the manual pages.
Time: 07:33:14 GMT, January 31, 2024