Polymake Template Library (PTL): pm Namespace Reference
Polymake Template Library (PTL)  4.2
pm Namespace Reference

global namespace for all classes from the polymake project More...

Namespaces

 AVL
 traits classes and such related to balanced trees
 

Classes

class  AccurateFloat
 minimalistic wrapper for MPFR numbers More...
 
class  Array
 Container class with constant time random access. More...
 
struct  attrib
 
class  Bitset
 Container class for dense sets of integers. More...
 
class  chunk_allocator
 Maintains a list of private memory chunks of fixed size. More...
 
class  color_error
 An exception of this type is thrown by an attempt to assign a wrong value to some color component. More...
 
class  Complement
 Complement as GenericSet. More...
 
class  conv
 Explicit type converter. More...
 
struct  Div
 result of integer division of two numbers (a,b) More...
 
class  EquivalenceRelation
 An equivalence relation on the integers 0,..,n-1 for a given size n. More...
 
struct  ExtGCD
 result of the extended gcd calculation for two numbers (a, b) More...
 
class  FaceMap
 
class  FacetList
 
struct  first_of_equal
 special tags for find_nearest denoting the first and last occurrence of a given key in a multi-set More...
 
struct  function_argument
 
class  GenericGraph
 Generic type for all graph classes. More...
 
class  GenericMatrix
 Generic type for matrices More...
 
class  GenericMutableSet
 Generic type for ordered mutable sets More...
 
class  GenericSet
 Generic type for ordered sets More...
 
class  GenericVector
 Generic type for vectors More...
 
class  Heap
 
class  HSV
 Color description in HSV space. More...
 
class  IncidenceMatrix
 0/1 incidence matrix. More...
 
class  input_truncator
 
class  Integer
 Integral number of unlimited precision. More...
 
struct  is_ordered
 check whether the default comparator operations::cmp can be used with the given type More...
 
struct  is_unordered_comparable
 check whether the comparator operations::cmp_unordered can be used with the given type More...
 
class  iterator_product
 
struct  list_search
 
struct  list_search_all
 
class  ListMatrix
 List of row vectors. More...
 
class  Map
 Associative array based on AVL::tree. More...
 
class  Matrix
 Matrix type class which holds the elements in a contiguous array
Additional arithmetic operations for matrices and useful constructions (unit_matrix, diag, ...) are listed at operations. More...
 
struct  merge_list
 
class  no_match
 
class  NormalRandom
 
struct  nothing
 Structure denoting the absence of data. More...
 
class  output_predicate_selector
 
class  permutation_iterator< permutations_heap >
 Implementation of the Heap's algorithm by R. Sedgewick. More...
 
class  PowerSet
 A Set with elements of type Set<E>, providing methods for adding elements while preserving subset or superset independence. Comparator is a functor defining a total ordering on Set<E>. More...
 
class  ptr_wrapper
 Wrapper for a pointer used as an iterator. More...
 
class  QuadraticExtension
 Realizes quadratic extensions of fields. More...
 
class  RandomPoints
 Generator of uniformly distributed random points on the unit sphere in R^d. More...
 
class  range_contractor
 
class  Rational
 Rational number with unlimited precision. More...
 
class  RGB
 Color description in RGB space: Red-Green-Blue additive color model. More...
 
class  Set
 An associative container based on a balanced binary search (AVL) tree. Comparator is a functor defining a total ordering on the element value domain. In most cases, the default choice (lexicographical order) will suffice for your needs. More...
 
class  Set_with_dim
 Set_with_dim as GenericSet. More...
 
class  shared_array
 
class  shared_object
 
class  SingleElementSetCmp
 A set consisting of exactly one element. More...
 
class  SmithNormalForm
 Complete result of computation of Smith normal form. More...
 
class  SparseMatrix
 A two-dimensional associative array with row and column indices as keys.
More...
 
class  SparseMatrixStatistics
 sparse matrix statistics collection More...
 
class  SparseVector
 
struct  spec_object_traits< TropicalNumber< Addition, Scalar > >
 
class  TropicalNumber
 
class  unary_predicate_selector
 
class  UniformlyRandom< AccurateFloat >
 Generator of random AccurateFloat numbers from [0, 1) More...
 
class  UniformlyRandom< Bitset >
 Generator of random Bitset of a given maximal cardinality. More...
 
class  UniformlyRandom< Rational >
 
struct  unlimited
 
class  Vector
 Vector type class which holds the elements in a contiguous array More...
 

Functions

Int incl (const Bitset &s1, const Bitset &s2) noexcept
 See incl(GenericSet,GenericSet)
 
template<typename Matrix1 , typename Matrix2 >
auto diag (const GenericIncidenceMatrix< Matrix1 > &m1, const GenericIncidenceMatrix< Matrix2 > &m2) -> decltype(make_block_diag< false >(m1.top(), m2.top()))
 create a block-diagonal incidence matrix
 
template<typename Matrix1 , typename Matrix2 >
auto diag_1 (const GenericIncidenceMatrix< Matrix1 > &m1, const GenericIncidenceMatrix< Matrix2 > &m2) -> decltype(make_block_diag< true >(m1.top(), m2.top()))
 create a block-diagonal incidence matrix, fill the corner blocks with 1's
 
template<typename TargetType , typename TMatrix >
const TMatrix & convert_to (const GenericMatrix< TMatrix, TargetType > &m)
 explicit conversion of matrix elements to another type
 
template<typename TVector , typename = std::enable_if_t<is_generic_vector<TVector>::value>>
auto repeat_row (TVector &&v, Int n=0) -> RepeatedRow< diligent_ref_t< unwary_t< TVector >>>
 Create a matrix with n rows, each equal to v.
 
template<typename TVector , typename = std::enable_if_t<is_generic_vector<TVector>::value>>
auto repeat_col (TVector &&v, Int n=0) -> RepeatedCol< diligent_ref_t< unwary_t< TVector >>>
 Create a matrix with n columns, each equal to v.
 
template<typename E >
auto same_element_matrix (E &&x, Int m, Int n)
 
template<typename E >
auto ones_matrix (Int m, Int n)
 
template<typename E >
auto zero_matrix (Int m, Int n)
 
template<typename TVector >
auto vector2row (GenericVector< TVector > &v)
 disguise a GenericVector as a matrix with 1 row
 
template<typename TVector >
auto vector2col (GenericVector< TVector > &v)
 disguise a GenericVector as a matrix with 1 column
 
template<typename TVector >
auto diag (const GenericVector< TVector > &v)
 Create a square diagonal matrix from a GenericVector.
 
template<typename TVector >
auto anti_diag (const GenericVector< TVector > &v)
 Create a anti-diagonal matrix.
 
template<typename E >
auto unit_matrix (Int dim)
 Create a unit_matrix of dimension dim.
 
template<typename E , typename Matrix1 , typename Matrix2 >
auto diag (const GenericMatrix< Matrix1, E > &m1, const GenericMatrix< Matrix2, E > &m2)
 Create a block-diagonal matrix.
 
template<typename E , typename Vector1 , typename Matrix2 >
auto diag (const GenericVector< Vector1, E > &v1, const GenericMatrix< Matrix2, E > &m2)
 
template<typename E , typename Matrix1 , typename Vector2 >
auto diag (const GenericMatrix< Matrix1, E > &m1, const GenericVector< Vector2, E > &v2)
 
template<typename E , typename Matrix1 , typename Matrix2 >
auto anti_diag (const GenericMatrix< Matrix1, E > &m1, const GenericMatrix< Matrix2, E > &m2)
 Create a block-anti-diagonal matrix.
 
template<typename E , typename Vector1 , typename Matrix2 >
auto anti_diag (const GenericVector< Vector1, E > &v1, const GenericMatrix< Matrix2, E > &m2)
 
template<typename E , typename Matrix1 , typename Vector2 >
auto anti_diag (const GenericMatrix< Matrix1, E > &m1, const GenericVector< Vector2, E > &v2)
 
template<typename Comparator , typename E >
auto scalar2set (E &&x)
 construct an one-element set with explicitly specified comparator type
 
template<typename E >
auto scalar2set (E &&x)
 construct an one-element set with standard (lexicographical) comparator type
 
template<typename Set1 , typename Set2 , typename E1 , typename E2 , class Comparator >
Int incl (const GenericSet< Set1, E1, Comparator > &s1, const GenericSet< Set2, E2, Comparator > &s2)
 
template<typename TargetType , typename TVector >
const TVector & convert_to (const GenericVector< TVector, TargetType > &v)
 explicit conversion of vector elements to another type
 
template<typename E >
auto same_element_vector (E &&x, Int dim=0)
 Create a vector with all entries equal to the given element x.
 
template<typename E >
auto ones_vector (Int dim=0)
 Create a vector with all entries equal to 1.
 
template<typename E >
auto zero_vector (Int dim=0)
 Create a vector with all entries equal to 0.
 
graph::Graph complete_graph (Int n_nodes)
 create complete graph on given number of nodes
 
template<typename Matrix1 , typename Matrix2 >
IncidenceMatrix convolute (const GenericIncidenceMatrix< Matrix1 > &m1, const GenericIncidenceMatrix< Matrix2 > &m2)
 Convolution of two incidence relations.
 
constexpr bool isfinite (const __mpz_struct &)
 data from third parties can't have infinite values
 
Int isinf (double x) noexcept
 
template<typename Target , typename Source >
inherit_reference_t< Target, Source && > convert_to (Source &&x, std::enable_if_t< is_derived_from< pure_type_t< Source >, Target >::value||std::is_same< pure_type_t< Source >, Target >::value, void ** >=nullptr)
 trivial conversion of an object to itself (mutable access)
 
template<typename Target , typename Source >
Target convert_to (const Source &x, std::enable_if_t< can_initialize< Source, Target >::value &&(!std::is_arithmetic< Source >::value||!std::is_arithmetic< Target >::value) &&!is_derived_from< Source, Target >::value, void ** >=nullptr)
 conversion via functor
 
template<typename E >
std::enable_if_t< is_field< E >::value, E > det (Matrix< E > M)
 determinant of a matrix
 
template<typename E >
std::enable_if_t< is_field< E >::value, Matrix< E > > inv (Matrix< E > M)
 matrix inversion
 
template<typename E >
std::enable_if_t< is_field< E >::value, Vector< E > > lin_solve (Matrix< E > A, Vector< E > b)
 solving systems of linear equations
 
template<typename Iterator , typename Predicate >
auto make_unary_predicate_selector (Iterator &&it, const Predicate &pred)
 Convenience function creating an iterator with element selection.
 
template<typename Matrix2 , typename Vector2 , typename E >
bool has_solution (const GenericMatrix< Matrix2, E > &A, const GenericVector< Vector2, E > &B)
 
template<typename T , typename... Args>
T * construct_at (T *place, Args &&... args)
 
template<typename Iterator >
void normalize (Iterator dst)
 Divide each vector in a sequence thru its length (L2-norm)
 
template<typename TMatrix , typename E >
Vector< E > barycenter (const GenericMatrix< TMatrix, E > &V)
 Compute the average over the rows of a matrix.
 
template<typename TMatrix , typename E >
std::enable_if_t< is_field< E >::value, E > det (const GenericMatrix< TMatrix, E > &m)
 Compute the determinant of a matrix using the Gauss elimination method.
 
template<typename TMatrix , typename E >
trace (const GenericMatrix< TMatrix, E > &m)
 Compute the trace of a matrix.
 
template<typename TMatrix , typename TVector , typename E >
std::enable_if_t< is_field< E >::value, typename TVector::persistent_type > reduce (const GenericMatrix< TMatrix, E > &A, const GenericVector< TVector, E > &V)
 Reduce a vector with a given matrix using the Gauss elimination method.
 
template<typename TMatrix , typename E >
std::enable_if_t< is_field< E >::value, typename TMatrix::persistent_nonsymmetric_type > inv (const GenericMatrix< TMatrix, E > &m)
 
template<typename TMatrix , typename TVector , typename E >
std::enable_if_t< is_field< E >::value, Vector< E > > lin_solve (const GenericMatrix< TMatrix, E > &A, const GenericVector< TVector, E > &b)
 
template<typename AHRowIterator , typename VectorType , typename RowBasisOutputIterator , typename ColBasisOutputIterator >
bool project_rest_along_row (AHRowIterator &h, const VectorType &v, RowBasisOutputIterator row_basis_consumer, ColBasisOutputIterator col_basis_consumer, Int i=0)
 
template<typename VectorType , typename RowBasisOutputIterator , typename ColBasisOutputIterator , typename E >
bool basis_of_rowspan_intersect_orthogonal_complement (ListMatrix< SparseVector< E > > &H, const VectorType &v, RowBasisOutputIterator row_basis_consumer, ColBasisOutputIterator col_basis_consumer, Int i=0)
 
template<typename T , typename R >
bool add_row_if_rowspace_increases (ListMatrix< SparseVector< T > > &M, const SparseVector< T > &v, ListMatrix< SparseVector< R > > &kernel_so_far)
 
template<typename VectorIterator , typename RowBasisOutputIterator , typename ColBasisOutputIterator , typename AH_matrix >
void null_space (VectorIterator v, RowBasisOutputIterator row_basis_consumer, ColBasisOutputIterator col_basis_consumer, AH_matrix &H, bool simplify=false)
 
template<typename TVector , typename E >
ListMatrix< SparseVector< E > > null_space_oriented (const GenericVector< TVector, E > &V, Int req_sign)
 
template<typename TMatrix , typename E >
std::pair< Set< Int >, Set< Int > > basis_affine (const GenericMatrix< TMatrix, E > &M)
 The same as basis(), but ignoring the first column of the matrix.
 
template<typename E , typename Vector1 , typename Vector2 >
Vector2::persistent_type proj (const GenericVector< Vector1, E > &u, const GenericVector< Vector2, E > &v)
 Project u to v.
 
template<typename Matrix1 , typename Matrix2 >
void project_to_orthogonal_complement (Matrix1 &M, const Matrix2 &N)
 
template<typename TVector >
Set< Int > support (const GenericVector< TVector > &v)
 the indices of nonzero entries
 
template<typename Vector1 , typename Vector2 >
Vector1::persistent_type reflect (const GenericVector< Vector1 > &u, const GenericVector< Vector2 > &nv)
 reflect u in the plane normal to nv
 
template<typename TVector >
TVector::persistent_type dehomogenize (const GenericVector< TVector > &V)
 
template<typename TMatrix >
TMatrix::persistent_nonsymmetric_type dehomogenize (const GenericMatrix< TMatrix > &M)
 
template<typename TVector >
TVector::persistent_type dehomogenize_trop (const GenericVector< TVector > &V)
 
template<typename TMatrix >
TMatrix::persistent_nonsymmetric_type dehomogenize_trop (const GenericMatrix< TMatrix > &M)
 
template<typename TMatrix >
TMatrix::persistent_nonsymmetric_type remove_zero_rows (const GenericMatrix< TMatrix > &m)
 Remove all matrix rows that contain only zeros.
 
template<typename VectorIterator >
void orthogonalize (VectorIterator v)
 Apply the Gram-Schmidt orthogonalization to the vector sequence.
 
template<typename VectorIterator >
void orthogonalize_affine (VectorIterator v)
 
template<typename TMatrix >
Set< Int > far_points (const GenericMatrix< TMatrix > &M)
 Find row indices of all far points (that is, having zero in the first column).
 
template<typename E , typename TMatrix , typename TVector >
Set< Int > orthogonal_rows (const GenericMatrix< TMatrix, E > &M, const GenericVector< TVector, E > &v)
 Find indices of rows orthogonal to the given vector.
 
template<typename Coefficient , typename Exponent >
UniPolynomial< Coefficient, Exponent > div_exact (const UniPolynomial< Coefficient, Exponent > &a, const UniPolynomial< Coefficient, Exponent > &b)
 
template<typename TargetType , typename Exponent >
const Polynomial< TargetType, Exponent > & convert_to (const Polynomial< TargetType, Exponent > &p)
 explicit conversion of polynomial coefficients to another type
 
template<typename Iterator >
PowerSet< typename iterator_traits< Iterator >::value_type::element_type, typename iterator_traits< Iterator >::value_type::element_comparator > ridges (Iterator set)
 Gather all independent intersections of subsets from the given PowerSet.
 
constexpr bool isfinite (const __mpq_struct &)
 data from third parties can't have infinite values
 
template<typename E >
Int smith_normal_form_only (SparseMatrix< E > &M, std::list< std::pair< E, Int >> &torsion)
 
template<typename Matrix , typename E >
SmithNormalForm< E > smith_normal_form (const GenericMatrix< Matrix, E > &M, std::enable_if_t< std::numeric_limits< E >::is_integer, bool > inverse_companions=false)
 

Detailed Description

global namespace for all classes from the polymake project

The most common classes/functions/... are exported to the namespace polymake. Use this for your client code.

Function Documentation

◆ add_row_if_rowspace_increases()

template<typename T , typename R >
bool pm::add_row_if_rowspace_increases ( ListMatrix< SparseVector< T > > &  M,
const SparseVector< T > &  v,
ListMatrix< SparseVector< R > > &  kernel_so_far 
)

add a row to a matrix M iff this increases the rowspan of M. As a side effect, update the kernel of M.

Parameters
Ma matrix, implemented as a ListMatrix so that the row addition is cheap
kernel_so_fara matrix whose rows are supposed to be orthogonal to the rows of M.
va vector to be added to M iff this increases the dimension of the rowspan of M. We allow the entries of v and kernel_so_far to be of different type, e.g. T=Integer, R=Rational

◆ basis_of_rowspan_intersect_orthogonal_complement()

template<typename VectorType , typename RowBasisOutputIterator , typename ColBasisOutputIterator , typename E >
bool pm::basis_of_rowspan_intersect_orthogonal_complement ( ListMatrix< SparseVector< E > > &  H,
const VectorType &  v,
RowBasisOutputIterator  row_basis_consumer,
ColBasisOutputIterator  col_basis_consumer,
Int  i = 0 
)

Compute a basis of (subspace spanned by the rows of H) intersect (orthogonal complement of v) Projects all row vectors in H onto the orthogonal complement of v, and keeps a basis of rowspan(H) intersect orthogonal(v)

Parameters
vinput iterator over the vectors
row_basis_consumeroutput iterator consuming the indices of basis rows (=vectors)
col_basis_consumeroutput iterator consuming the indices of basis columns (=coordinates)
Returns
modified_H boolean value indicating whether the kernel matrix H has been modified

◆ construct_at()

template<typename T , typename... Args>
T* pm::construct_at ( T *  place,
Args &&...  args 
)

These functions are defined in the standard libraries at non-portable locations under non-portable names. Until they are unified, it's easier to provide an own implementation.

◆ dehomogenize() [1/2]

template<typename TMatrix >
TMatrix::persistent_nonsymmetric_type pm::dehomogenize ( const GenericMatrix< TMatrix > &  M)

Divide rowwise by the elements of the first column and strip it off. As a special case, an empty matrix is passed unchanged.

◆ dehomogenize() [2/2]

template<typename TVector >
TVector::persistent_type pm::dehomogenize ( const GenericVector< TVector > &  V)

Divide by the first element and strip it off. As a special case, an empty vector is passed unchanged.

◆ dehomogenize_trop() [1/2]

template<typename TMatrix >
TMatrix::persistent_nonsymmetric_type pm::dehomogenize_trop ( const GenericMatrix< TMatrix > &  M)

Subtract rowwise the elements of the first column and strip it off. As a special case, an empty matrix is passed unchanged.

◆ dehomogenize_trop() [2/2]

template<typename TVector >
TVector::persistent_type pm::dehomogenize_trop ( const GenericVector< TVector > &  V)

Subtract the first element and strip it off. As a special case, an empty vector is passed unchanged.

◆ div_exact()

template<typename Coefficient , typename Exponent >
UniPolynomial<Coefficient, Exponent> pm::div_exact ( const UniPolynomial< Coefficient, Exponent > &  a,
const UniPolynomial< Coefficient, Exponent > &  b 
)

Perform the polynomial division, discarding the remainder. Although the name suggests that the divisor must be a factor of the divident, you can put arbitrary polynomials here. The name is rather chosen for compatibility with the Integer class.

◆ has_solution()

template<typename Matrix2 , typename Vector2 , typename E >
bool pm::has_solution ( const GenericMatrix< Matrix2, E > &  A,
const GenericVector< Vector2, E > &  B 
)
inline

Solve the linear system A*x==B

Returns
x
Exceptions
degenerate_matrixif det(A) == 0
infeasibleif rank(A) != rank(A|B)

◆ inv()

template<typename TMatrix , typename E >
std::enable_if_t<is_field<E>::value, typename TMatrix::persistent_nonsymmetric_type> pm::inv ( const GenericMatrix< TMatrix, E > &  m)

Compute the inverse matrix $A^-1$ using the Gauss elimination method.

Exceptions
degenerate_matrixif det(A)==0

◆ isinf()

Int pm::isinf ( double  x)
inlinenoexcept

return the sign of the inifinite value, or 0 if the value is finite std::isinf returns bool nowadays which is insufficient for efficient operations

◆ lin_solve()

template<typename TMatrix , typename TVector , typename E >
std::enable_if_t<is_field<E>::value, Vector<E> > pm::lin_solve ( const GenericMatrix< TMatrix, E > &  A,
const GenericVector< TVector, E > &  b 
)

Solve the linear system A*x==b

Returns
xapps/polytope/src/reverse_search_graph.cco
Exceptions
degenerate_matrixif det(A) == 0
infeasibleif rank(A) != rank(A|b)

◆ null_space()

template<typename VectorIterator , typename RowBasisOutputIterator , typename ColBasisOutputIterator , typename AH_matrix >
void pm::null_space ( VectorIterator  v,
RowBasisOutputIterator  row_basis_consumer,
ColBasisOutputIterator  col_basis_consumer,
AH_matrix &  H,
bool  simplify = false 
)

Compute the basis of the subspaces spanned by a sequence of vectors and orthogonal one. Projects all row vectors in H into the orthogonal complement of the vectors that v iterates over, and keeps a basis of rowspan(H) intersect orthogonal(v_i)

Parameters
vinput iterator over the vectors
row_basis_consumeroutput iterator consuming the indices of basis rows (=vectors)
col_basis_consumeroutput iterator consuming the indices of basis columns (=coordinates)

◆ null_space_oriented()

template<typename TVector , typename E >
ListMatrix< SparseVector<E> > pm::null_space_oriented ( const GenericVector< TVector, E > &  V,
Int  req_sign 
)
Parameters
req_signexpected sign of det( null_space(V) / V )

◆ orthogonalize_affine()

template<typename VectorIterator >
void pm::orthogonalize_affine ( VectorIterator  v)

The same as orthogonalize(.), but making the affine parts of the resulting vectors (without 0-th coordinate) orthogonal

◆ project_rest_along_row()

template<typename AHRowIterator , typename VectorType , typename RowBasisOutputIterator , typename ColBasisOutputIterator >
bool pm::project_rest_along_row ( AHRowIterator &  h,
const VectorType &  v,
RowBasisOutputIterator  row_basis_consumer,
ColBasisOutputIterator  col_basis_consumer,
Int  i = 0 
)

if the row referenced by h is not orthogonal to v, project out *h from all subsequent rows so that they become orthogonal to v

Parameters
hrow iterator over a matrix
vthe vector that the rows should become orthogonal to
row_basis_consumeroutput iterator consuming the indices of basis rows (=vectors)
col_basis_consumeroutput iterator consuming the indices of basis columns (=coordinates)
Returns
true if the matrix has been modified, false otherwise

◆ project_to_orthogonal_complement()

template<typename Matrix1 , typename Matrix2 >
void pm::project_to_orthogonal_complement ( Matrix1 &  M,
const Matrix2 &  N 
)

project the rows of M into the orthogonal complement of N the rows of N need to be orthogonal

◆ smith_normal_form()

template<typename Matrix , typename E >
SmithNormalForm<E> pm::smith_normal_form ( const GenericMatrix< Matrix, E > &  M,
std::enable_if_t< std::numeric_limits< E >::is_integer, bool >  inverse_companions = false 
)

Compute the Smith normal form and companion matrices.

Parameters
inverse_companionsif true: result.form = result.left_companion * M * result.right_companion if false: M = result.left_companion * result.form * result.right_companion

◆ smith_normal_form_only()

template<typename E >
Int pm::smith_normal_form_only ( SparseMatrix< E > &  M,
std::list< std::pair< E, Int >> &  torsion 
)

Compute the compact respresentation of the Smith normal form, without companion matrices. Input matrix M is corrupted during the computations.

Parameters
[out]torsionlist of diagonal elements of SNF not equal 1 with multiplicities.
Returns
rank of M.