Title: | Volume Approximation and Sampling of Convex Polytopes |
---|---|
Description: | Provides an R interface for 'volesti' C++ package. 'volesti' computes estimations of volume of polytopes given by (i) a set of points, (ii) linear inequalities or (iii) Minkowski sum of segments (a.k.a. zonotopes). There are three algorithms for volume estimation as well as algorithms for sampling, rounding and rotating polytopes. Moreover, 'volesti' provides algorithms for estimating copulas useful in computational finance. Methods implemented in 'volesti' are described in A. Chalkis and V. Fisikopoulos (2022) <doi:10.32614/RJ-2021-077> and references therein. |
Authors: | Vissarion Fisikopoulos [aut, cre, cph] , Apostolos Chalkis [aut, cph] |
Maintainer: | Vissarion Fisikopoulos <[email protected]> |
License: | LGPL-3 |
Version: | 1.1.2-8 |
Built: | 2024-11-09 03:52:18 UTC |
Source: | https://github.com/cran/volesti |
Given a matrix that contains row-wise the assets' returns and a sliding window win_length
, this function computes an approximation of the joint distribution (copula, e.g. see https://en.wikipedia.org/wiki/Copula_(probability_theory)) between portfolios' return and volatility in each time period defined by win_len
.
For each copula it computes an indicator: If the indicator is large it corresponds to a crisis period and if it is small it corresponds to a normal period.
In particular, the periods over which the indicator is greater than 1 for more than 60 consecutive sliding windows are warnings and for more than 100 are crisis. The sliding window is shifted by one day.
compute_indicators( returns, parameters = list(win_length = 60, m = 100, n = 5e+05, nwarning = 60, ncrisis = 100) )
compute_indicators( returns, parameters = list(win_length = 60, m = 100, n = 5e+05, nwarning = 60, ncrisis = 100) )
returns |
A |
parameters |
A list to set a parameterization.
|
A list that contains the indicators and the corresponding vector that label each time period with respect to the market state: a) normal, b) crisis, c) warning.
L. Cales, A. Chalkis, I.Z. Emiris, V. Fisikopoulos, “Practical volume computation of structured convex bodies, and an application to modeling portfolio dependencies and financial crises,” Proc. of Symposium on Computational Geometry, Budapest, Hungary, 2018.
# simple example on random asset returns asset_returns = replicate(10, rnorm(14)) market_states_and_indicators = compute_indicators(asset_returns, parameters = list("win_length" = 10, "m" = 10, "n" = 10000, "nwarning" = 2, "ncrisis" = 3))
# simple example on random asset returns asset_returns = replicate(10, rnorm(14)) market_states_and_indicators = compute_indicators(asset_returns, parameters = list("win_length" = 10, "m" = 10, "n" = 10000, "nwarning" = 2, "ncrisis" = 3))
Given two families of parallel hyperplanes or a family of parallel hyperplanes and a family of concentric ellispoids centered at the origin intersecting the canonical simplex, this function uniformly samples from the canonical simplex and construct an approximation of the bivariate probability distribution, called copula (see https://en.wikipedia.org/wiki/Copula_(probability_theory)). At least two families of hyperplanes or one family of hyperplanes and one family of ellipsoids have to be given as input.
copula(r1, r2 = NULL, sigma = NULL, m = NULL, n = NULL, seed = NULL)
copula(r1, r2 = NULL, sigma = NULL, m = NULL, n = NULL, seed = NULL)
r1 |
The |
r2 |
Optional. The |
sigma |
Optional. The |
m |
The number of the slices for the copula. The default value is 100. |
n |
The number of points to sample. The default value is |
seed |
Optional. A fixed seed for the number generator. |
A numerical matrix that corresponds to a copula.
L. Cales, A. Chalkis, I.Z. Emiris, V. Fisikopoulos, “Practical volume computation of structured convex bodies, and an application to modeling portfolio dependencies and financial crises,” Proc. of Symposium on Computational Geometry, Budapest, Hungary, 2018.
# compute a copula for two random families of parallel hyperplanes h1 = runif(n = 10, min = 1, max = 1000) h1 = h1 / 1000 h2=runif(n = 10, min = 1, max = 1000) h2 = h2 / 1000 cop = copula(r1 = h1, r2 = h2, m = 10, n = 100000) # compute a copula for a family of parallel hyperplanes and a family of conentric ellipsoids h = runif(n = 10, min = 1, max = 1000) h = h / 1000 E = replicate(10, rnorm(20)) E = cov(E) cop = copula(r1 = h, sigma = E, m = 10, n = 100000)
# compute a copula for two random families of parallel hyperplanes h1 = runif(n = 10, min = 1, max = 1000) h1 = h1 / 1000 h2=runif(n = 10, min = 1, max = 1000) h2 = h2 / 1000 cop = copula(r1 = h1, r2 = h2, m = 10, n = 100000) # compute a copula for a family of parallel hyperplanes and a family of conentric ellipsoids h = runif(n = 10, min = 1, max = 1000) h = h / 1000 E = replicate(10, rnorm(20)) E = cov(E) cop = copula(r1 = h, sigma = E, m = 10, n = 100000)
The -dimensional unit simplex is the set of points
, s.t.:
,
. The
-dimensional canonical simplex is the set of points
, s.t.:
,
.
direct_sampling(body, n)
direct_sampling(body, n)
body |
A list to request exact uniform sampling from special well known convex bodies through the following input parameters:
|
n |
The number of points that the function is going to sample. |
A matrix that contains, column-wise, the sampled points from the convex polytope P.
R.Y. Rubinstein and B. Melamed, “Modern simulation and modeling” Wiley Series in Probability and Statistics, 1998.
A Smith, Noah and W Tromble, Roy, “Sampling Uniformly from the Unit Simplex,” Center for Language and Speech Processing Johns Hopkins University, 2004.
# 100 uniform points from the 2-d unit ball points = direct_sampling(n = 100, body = list("type" = "ball", "dimension" = 2))
# 100 uniform points from the 2-d unit ball points = direct_sampling(n = 100, body = list("type" = "ball", "dimension" = 2))
Given a zonotope (as an object of class Zonotope), this function computes the sum of the absolute values of the determinants of all the submatrices of the
matrix
that contains row-wise the
-dimensional segments that define the zonotope.
For an arbitrary simplex that is given in V-representation this function computes the absolute value of the determinant formed by the simplex's points assuming it is shifted to the origin.
exact_vol(P)
exact_vol(P)
P |
A polytope |
The exact volume of the input polytope, for zonotopes, simplices in V-representation and polytopes with known exact volume
E. Gover and N. Krikorian, “Determinants and the Volumes of Parallelotopes and Zonotopes,” Linear Algebra and its Applications, 433(1), 28 - 40, 2010.
# compute the exact volume of a 5-dimensional zonotope defined by the Minkowski sum of 10 segments Z = gen_rand_zonotope(2, 5) vol = exact_vol(Z) # compute the exact volume of a 2-d arbitrary simplex V = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Vpolytope(V = V) vol = exact_vol(P) # compute the exact volume the 10-dimensional cross polytope P = gen_cross(10,'V') vol = exact_vol(P)
# compute the exact volume of a 5-dimensional zonotope defined by the Minkowski sum of 10 segments Z = gen_rand_zonotope(2, 5) vol = exact_vol(Z) # compute the exact volume of a 2-d arbitrary simplex V = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Vpolytope(V = V) vol = exact_vol(P) # compute the exact volume the 10-dimensional cross polytope P = gen_cross(10,'V') vol = exact_vol(P)
A half-space is given as a pair of a vector
and a scalar
s.t.:
. This function calls the Ali's version of the Varsi formula to compute a frustum of the simplex.
frustum_of_simplex(a, z0)
frustum_of_simplex(a, z0)
a |
A |
z0 |
The scalar that defines the half-space. |
The percentage of the volume of the simplex that is contained in the intersection of a given half-space and the simplex.
Varsi, Giulio, “The multidimensional content of the frustum of the simplex,” Pacific J. Math. 46, no. 1, 303–314, 1973.
Ali, Mir M., “Content of the frustum of a simplex,” Pacific J. Math. 48, no. 2, 313–322, 1973.
# compute the frustum of H: -x1+x2<=0 a=c(-1,1) z0=0 frustum = frustum_of_simplex(a, z0)
# compute the frustum of H: -x1+x2<=0 a=c(-1,1) z0=0 frustum = frustum_of_simplex(a, z0)
This function generates the -dimensional cross polytope in H- or V-representation.
gen_cross(dimension, representation = "H")
gen_cross(dimension, representation = "H")
dimension |
The dimension of the cross polytope. |
representation |
A string to declare the representation. It has to be |
A polytope class representing a cross polytope in H- or V-representation.
# generate a 10-dimensional cross polytope in H-representation P = gen_cross(5, 'H') # generate a 15-dimension cross polytope in V-representation P = gen_cross(15, 'V')
# generate a 10-dimensional cross polytope in H-representation P = gen_cross(5, 'H') # generate a 15-dimension cross polytope in V-representation P = gen_cross(15, 'V')
This function generates the -dimensional unit hypercube
in H- or V-representation.
gen_cube(dimension, representation = "H")
gen_cube(dimension, representation = "H")
dimension |
The dimension of the hypercube |
representation |
A string to declare the representation. It has to be |
A polytope class representing the unit -dimensional hypercube in H- or V-representation.
# generate a 10-dimensional hypercube in H-representation P = gen_cube(10, 'H') # generate a 15-dimension hypercube in V-representation P = gen_cube(5, 'V')
# generate a 10-dimensional hypercube in H-representation P = gen_cube(10, 'H') # generate a 15-dimension hypercube in V-representation P = gen_cube(5, 'V')
This function generates a -dimensional polytope that is defined as the product of two
-dimensional unit simplices in H-representation.
gen_prod_simplex(dimension)
gen_prod_simplex(dimension)
dimension |
The dimension of the simplices. |
A polytope class representing the product of the two -dimensional unit simplices in H-representation.
# generate a product of two 5-dimensional simplices. P = gen_prod_simplex(5)
# generate a product of two 5-dimensional simplices. P = gen_prod_simplex(5)
This function generates a -dimensional polytope in H-representation with
facets. We pick
random hyperplanes tangent on the
-dimensional unit hypersphere as facets.
gen_rand_hpoly(dimension, nfacets, generator = list(constants = "sphere"))
gen_rand_hpoly(dimension, nfacets, generator = list(constants = "sphere"))
dimension |
The dimension of the convex polytope. |
nfacets |
The number of the facets. |
generator |
A list that could contain two elements.
|
A polytope class representing a H-polytope.
# generate a 10-dimensional polytope with 50 facets P = gen_rand_hpoly(10, 50)
# generate a 10-dimensional polytope with 50 facets P = gen_rand_hpoly(10, 50)
This function generates a -dimensional polytope in V-representation with
vertices. We pick
random points from the boundary of the
-dimensional unit hypersphere as vertices.
gen_rand_vpoly(dimension, nvertices, generator = list(body = "sphere"))
gen_rand_vpoly(dimension, nvertices, generator = list(body = "sphere"))
dimension |
The dimension of the convex polytope. |
nvertices |
The number of the vertices. |
generator |
A list that could contain two elements.
|
A polytope class representing a V-polytope.
# generate a 10-dimensional polytope defined as the convex hull of 25 random vertices P = gen_rand_vpoly(10, 25)
# generate a 10-dimensional polytope defined as the convex hull of 25 random vertices P = gen_rand_vpoly(10, 25)
This function generates a random -dimensional zonotope defined by the Minkowski sum of
-dimensional segments.
The function considers
random directions in
. There are three strategies to pick the length of each segment: a) it is uniformly sampled from
, b) it is random from
truncated to
, c) it is random from
truncated to
.
gen_rand_zonotope( dimension, nsegments, generator = list(distribution = "uniform") )
gen_rand_zonotope( dimension, nsegments, generator = list(distribution = "uniform") )
dimension |
The dimension of the zonotope. |
nsegments |
The number of segments that generate the zonotope. |
generator |
A list that could contain two elements.
|
A polytope class representing a zonotope.
# generate a 10-dimensional zonotope defined by the Minkowski sum of 20 segments P = gen_rand_zonotope(10, 20)
# generate a 10-dimensional zonotope defined by the Minkowski sum of 20 segments P = gen_rand_zonotope(10, 20)
This function generates the -dimensional unit simplex in H- or V-representation.
gen_simplex(dimension, representation = "H")
gen_simplex(dimension, representation = "H")
dimension |
The dimension of the unit simplex. |
representation |
A string to declare the representation. It has to be |
A polytope class representing the -dimensional unit simplex in H- or V-representation.
# generate a 10-dimensional simplex in H-representation PolyList = gen_simplex(10, 'H') # generate a 20-dimensional simplex in V-representation P = gen_simplex(20, 'V')
# generate a 10-dimensional simplex in H-representation PolyList = gen_simplex(10, 'H') # generate a 20-dimensional simplex in V-representation P = gen_simplex(20, 'V')
This function generates a -dimensional skinny hypercube
.
gen_skinny_cube(dimension)
gen_skinny_cube(dimension)
dimension |
The dimension of the skinny hypercube. |
A polytope class representing the -dimensional skinny hypercube in H-representation.
# generate a 10-dimensional skinny hypercube. P = gen_skinny_cube(10)
# generate a 10-dimensional skinny hypercube. P = gen_skinny_cube(10)
An H-polytope is a convex polytope defined by a set of linear inequalities or equivalently a -dimensional H-polytope with
facets is defined by a
matrix A and a
-dimensional vector b, s.t.:
.
An numerical matrix.
An -dimensional vector b.
The volume of the polytope if it is known, otherwise by default.
A character with default value 'Hpolytope', to declare the representation of the polytope.
A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b)
A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b)
For a H-polytope described by a matrix
and a
-dimensional vector
, s.t.:
, this function computes the largest inscribed ball (Chebychev ball) by solving the corresponding linear program.
For both zonotopes and V-polytopes the function computes the minimum
s.t.:
for all
. Then the ball centered at the origin with radius
is an inscribed ball.
inner_ball(P)
inner_ball(P)
P |
A convex polytope. It is an object from class (a) Hpolytope or (b) Vpolytope or (c) Zonotope or (d) VpolytopeIntersection. |
A -dimensional vector that describes the inscribed ball. The first
coordinates corresponds to the center of the ball and the last one to the radius.
# compute the Chebychev ball of the 2d unit simplex P = gen_simplex(2,'H') ball_vec = inner_ball(P) # compute an inscribed ball of the 3-dimensional unit cube in V-representation P = gen_cube(3, 'V') ball_vec = inner_ball(P)
# compute the Chebychev ball of the 2d unit simplex P = gen_simplex(2,'H') ball_vec = inner_ball(P) # compute an inscribed ball of the 3-dimensional unit cube in V-representation P = gen_cube(3, 'V') ball_vec = inner_ball(P)
Read a SDPA format file and return a spectrahedron (an object of class Spectrahedron) which is defined by the linear matrix inequality in the input file, and the objective function.
read_sdpa_format_file(path)
read_sdpa_format_file(path)
path |
Name of the input file |
A list with two named items: an item "matrices" which is an object of class Spectrahedron and an vector "objFunction"
path = system.file('extdata', package = 'volesti') l = read_sdpa_format_file(paste0(path,'/sdpa_n2m3.txt')) Spectrahedron = l$spectrahedron objFunction = l$objFunction
path = system.file('extdata', package = 'volesti') l = read_sdpa_format_file(paste0(path,'/sdpa_n2m3.txt')) Spectrahedron = l$spectrahedron objFunction = l$objFunction
Given a convex H- or V- polytope or a zonotope or an intersection of two V-polytopes as input, this function applies (a) a random rotation or (b) a given rotation by an input matrix .
rotate_polytope(P, rotation = list())
rotate_polytope(P, rotation = list())
P |
A convex polytope. It is an object from class (a) Hpolytope, (b) Vpolytope, (c) Zonotope, (d) intersection of two V-polytopes. |
rotation |
A list that contains (a) the rotation matrix T and (b) the 'seed' to set a spesific seed for the number generator. |
Let be the given polytope and
the rotated one and
be the matrix of the linear transformation.
If is in H-representation and
is the matrix that contains the normal vectors of the facets of
then
contains the normal vactors of the facets of
.
If is in V-representation and
is the matrix that contains column-wise the vertices of
then
contains the vertices of
.
If is a zonotope and
is the matrix that contains column-wise the generators of
then
contains the generators of
.
If is a matrix that contains column-wise points in
then
contains points in
.
A list that contains the rotated polytope and the matrix of the linear transformation.
# rotate a H-polytope (2d unit simplex) P = gen_simplex(2, 'H') poly_matrix_list = rotate_polytope(P) # rotate a V-polytope (3d cube) P = gen_cube(3, 'V') poly_matrix_list = rotate_polytope(P) # rotate a 5-dimensional zonotope defined by the Minkowski sum of 15 segments Z = gen_rand_zonotope(3, 6) poly_matrix_list = rotate_polytope(Z)
# rotate a H-polytope (2d unit simplex) P = gen_simplex(2, 'H') poly_matrix_list = rotate_polytope(P) # rotate a V-polytope (3d cube) P = gen_cube(3, 'V') poly_matrix_list = rotate_polytope(P) # rotate a 5-dimensional zonotope defined by the Minkowski sum of 15 segments Z = gen_rand_zonotope(3, 6) poly_matrix_list = rotate_polytope(Z)
Given a convex H or V polytope or a zonotope as input this function brings the polytope in rounded position based on minimum volume enclosing ellipsoid of a pointset.
round_polytope(P, settings = list())
round_polytope(P, settings = list())
P |
A convex polytope. It is an object from class (a) Hpolytope or (b) Vpolytope or (c) Zonotope. |
settings |
Optional. A list to parameterize the method by the random walk.
|
A list with 4 elements: (a) a polytope of the same class as the input polytope class and (b) the element "T" which is the matrix of the inverse linear transformation that is applied on the input polytope, (c) the element "shift" which is the opposite vector of that which has shifted the input polytope, (d) the element "round_value" which is the determinant of the square matrix of the linear transformation that is applied on the input polytope.
I.Z.Emiris and V. Fisikopoulos, “Practical polytope volume approximation,” ACM Trans. Math. Soft., 2018.,
Michael J. Todd and E. Alper Yildirim, “On Khachiyan’s Algorithm for the Computation of Minimum Volume Enclosing Ellipsoids,” Discrete Applied Mathematics, 2007.
# rotate a H-polytope (2d unit simplex) A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b) listHpoly = round_polytope(P) # rotate a V-polytope (3d unit cube) using Random Directions HnR with step equal to 50 P = gen_cube(3, 'V') ListVpoly = round_polytope(P) # round a 2-dimensional zonotope defined by 6 generators using ball walk Z = gen_rand_zonotope(2,6) ListZono = round_polytope(Z)
# rotate a H-polytope (2d unit simplex) A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b) listHpoly = round_polytope(P) # rotate a V-polytope (3d unit cube) using Random Directions HnR with step equal to 50 P = gen_cube(3, 'V') ListVpoly = round_polytope(P) # round a 2-dimensional zonotope defined by 6 generators using ball walk Z = gen_rand_zonotope(2,6) ListZono = round_polytope(Z)
Sample n points with uniform or multidimensional spherical gaussian -with a mode at any point- as the target distribution.
sample_points(P, n, random_walk = NULL, distribution = NULL)
sample_points(P, n, random_walk = NULL, distribution = NULL)
P |
A convex polytope. It is an object from class (a) Hpolytope or (b) Vpolytope or (c) Zonotope or (d) VpolytopeIntersection. |
n |
The number of points that the function is going to sample from the convex polytope. |
random_walk |
Optional. A list that declares the random walk and some related parameters as follows:
|
distribution |
Optional. A list that declares the target density and some related parameters as follows:
|
A matrix that contains, column-wise, the sampled points from the convex polytope P.
# uniform distribution from the 3d unit cube in H-representation using ball walk P = gen_cube(3, 'H') points = sample_points(P, n = 100, random_walk = list("walk" = "BaW", "walk_length" = 5)) # gaussian distribution from the 2d unit simplex in H-representation with variance = 2 A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b) points = sample_points(P, n = 100, distribution = list("density" = "gaussian", "variance" = 2)) # uniform points from the boundary of a 2-dimensional random H-polytope P = gen_rand_hpoly(2,20) points = sample_points(P, n = 100, random_walk = list("walk" = "BRDHR"))
# uniform distribution from the 3d unit cube in H-representation using ball walk P = gen_cube(3, 'H') points = sample_points(P, n = 100, random_walk = list("walk" = "BaW", "walk_length" = 5)) # gaussian distribution from the 2d unit simplex in H-representation with variance = 2 A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE) b = c(0,0,1) P = Hpolytope(A = A, b = b) points = sample_points(P, n = 100, distribution = list("density" = "gaussian", "variance" = 2)) # uniform points from the boundary of a 2-dimensional random H-polytope P = gen_rand_hpoly(2,20) points = sample_points(P, n = 100, random_walk = list("walk" = "BRDHR"))
A spectrahedron is a convex body defined by a linear matrix inequality of the form .
The matrices
are symmetric
real matrices and
denoted negative semidefiniteness.
A List that contains the matrices .
A0 = matrix(c(-1,0,0,0,-2,1,0,1,-2), nrow=3, ncol=3, byrow = TRUE) A1 = matrix(c(-1,0,0,0,0,1,0,1,0), nrow=3, ncol=3, byrow = TRUE) A2 = matrix(c(0,0,-1,0,0,0,-1,0,0), nrow=3, ncol=3, byrow = TRUE) lmi = list(A0, A1, A2) S = Spectrahedron(matrices = lmi);
A0 = matrix(c(-1,0,0,0,-2,1,0,1,-2), nrow=3, ncol=3, byrow = TRUE) A1 = matrix(c(-1,0,0,0,0,1,0,1,0), nrow=3, ncol=3, byrow = TRUE) A2 = matrix(c(0,0,-1,0,0,0,-1,0,0), nrow=3, ncol=3, byrow = TRUE) lmi = list(A0, A1, A2) S = Spectrahedron(matrices = lmi);
For the volume approximation can be used three algorithms. Either CoolingBodies (CB) or SequenceOfBalls (SOB) or CoolingGaussian (CG). An H-polytope with facets is described by a
matrix
and a
-dimensional vector
, s.t.:
. A V-polytope is defined as the convex hull of
-dimensional points which correspond to the vertices of P. A zonotope is desrcibed by the Minkowski sum of
-dimensional segments.
volume(P, settings = NULL, rounding = FALSE)
volume(P, settings = NULL, rounding = FALSE)
P |
A convex polytope. It is an object from class a) Hpolytope or b) Vpolytope or c) Zonotope or d) VpolytopeIntersection. |
settings |
Optional. A list that declares which algorithm, random walk and values of parameters to use, as follows:
|
rounding |
A boolean parameter for rounding. The default value is |
The approximation of the volume of a convex polytope.
I.Z.Emiris and V. Fisikopoulos, “Practical polytope volume approximation,” ACM Trans. Math. Soft., 2018.,
A. Chalkis and I.Z.Emiris and V. Fisikopoulos, “Practical Volume Estimation by a New Annealing Schedule for Cooling Convex Bodies,” CoRR, abs/1905.05494, 2019.,
B. Cousins and S. Vempala, “A practical volume algorithm,” Springer-Verlag Berlin Heidelberg and The Mathematical Programming Society, 2015.
# calling SOB algorithm for a H-polytope (3d unit simplex) HP = gen_cube(3,'H') vol = volume(HP) # calling CG algorithm for a V-polytope (2d simplex) VP = gen_simplex(2,'V') vol = volume(VP, settings = list("algorithm" = "CG")) # calling CG algorithm for a 2-dimensional zonotope defined as the Minkowski sum of 4 segments Z = gen_rand_zonotope(2, 4) vol = volume(Z, settings = list("random_walk" = "RDHR", "walk_length" = 2))
# calling SOB algorithm for a H-polytope (3d unit simplex) HP = gen_cube(3,'H') vol = volume(HP) # calling CG algorithm for a V-polytope (2d simplex) VP = gen_simplex(2,'V') vol = volume(VP, settings = list("algorithm" = "CG")) # calling CG algorithm for a 2-dimensional zonotope defined as the Minkowski sum of 4 segments Z = gen_rand_zonotope(2, 4) vol = volume(Z, settings = list("random_walk" = "RDHR", "walk_length" = 2))
A V-polytope is a convex polytope defined by the set of its vertices.
An numerical matrix that contains the vertices row-wise.
The volume of the polytope if it is known, otherwise by default.
A character with default value 'Vpolytope', to declare the representation of the polytope.
V = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Vpolytope(V = V)
V = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Vpolytope(V = V)
An intersection of two V-polytopes is defined by the intersection of the two coresponding convex hulls.
An numerical matrix that contains the vertices of the first V-polytope (row-wise).
An numerical matrix that contains the vertices of the second V-polytope (row-wise).
The volume of the polytope if it is known, otherwise by default.
A character with default value 'VpolytopeIntersection', to declare the representation of the polytope.
P1 = gen_simplex(2,'V') P2 = gen_cross(2,'V') P = VpolytopeIntersection(V1 = P1@V, V2 = P2@V)
P1 = gen_simplex(2,'V') P2 = gen_cross(2,'V') P = VpolytopeIntersection(V1 = P1@V, V2 = P2@V)
Outputs a spectrahedron (the matrices defining a linear matrix inequality) and a vector (the objective function) to a SDPA format file.
write_sdpa_format_file(spectrahedron, objective_function, output_file)
write_sdpa_format_file(spectrahedron, objective_function, output_file)
spectrahedron |
A spectrahedron in n dimensions; must be an object of class Spectrahedron |
objective_function |
A numerical vector of length n |
output_file |
Name of the output file |
A0 = matrix(c(-1,0,0,0,-2,1,0,1,-2), nrow=3, ncol=3, byrow = TRUE) A1 = matrix(c(-1,0,0,0,0,1,0,1,0), nrow=3, ncol=3, byrow = TRUE) A2 = matrix(c(0,0,-1,0,0,0,-1,0,0), nrow=3, ncol=3, byrow = TRUE) lmi = list(A0, A1, A2) S = Spectrahedron(matrices = lmi) objFunction = c(1,1) write_sdpa_format_file(S, objFunction, "output.txt")
A0 = matrix(c(-1,0,0,0,-2,1,0,1,-2), nrow=3, ncol=3, byrow = TRUE) A1 = matrix(c(-1,0,0,0,0,1,0,1,0), nrow=3, ncol=3, byrow = TRUE) A2 = matrix(c(0,0,-1,0,0,0,-1,0,0), nrow=3, ncol=3, byrow = TRUE) lmi = list(A0, A1, A2) S = Spectrahedron(matrices = lmi) objFunction = c(1,1) write_sdpa_format_file(S, objFunction, "output.txt")
For the evaluation of the PCA method the exact volume of the approximation body is computed and the volume of the input zonotope is computed by CoolingBodies algorithm. The ratio of fitness is , where
is the approximate polytope.
zonotope_approximation( Z, fit_ratio = FALSE, settings = list(error = 0.1, walk_length = 1, win_len = 250, hpoly = FALSE) )
zonotope_approximation( Z, fit_ratio = FALSE, settings = list(error = 0.1, walk_length = 1, win_len = 250, hpoly = FALSE) )
Z |
A zonotope. |
fit_ratio |
Optional. A boolean parameter to request the computation of the ratio of fitness. |
settings |
Optional. A list that declares the values of the parameters of CB algorithm as follows:
|
A list that contains the approximation body in H-representation and the ratio of fitness
A.K. Kopetzki and B. Schurmann and M. Althoff, “Methods for Order Reduction of Zonotopes,” IEEE Conference on Decision and Control, 2017.
# over-approximate a 2-dimensional zonotope with 10 generators and compute the ratio of fitness Z = gen_rand_zonotope(2, 8) retList = zonotope_approximation(Z = Z)
# over-approximate a 2-dimensional zonotope with 10 generators and compute the ratio of fitness Z = gen_rand_zonotope(2, 8) retList = zonotope_approximation(Z = Z)
A zonotope is a convex polytope defined by the Minkowski sum of
-dimensional segments.
An numerical matrix that contains the segments (or generators) row-wise
The volume of the polytope if it is known, otherwise by default.
A character with default value 'Zonotope', to declare the representation of the polytope.
G = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Zonotope(G = G)
G = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE) P = Zonotope(G = G)