|
CGAL 6.1 - Homological Discrete Vector Fields
|
#include <CGAL/HDVF/Hdvf_core.h>
The class Hdvf_core is the core implementation of homological discrete vector fields (HDVF for short).
An enumeration FlagType is defined in the HDVF namespace and the Hdvf_core class maps each cell to one of the flags (namely PRIMARY, SECONDARY, CRITICAL). The NONE flag is used in child classes (such as Hdvf_duality or Hdvf_reduced) when computing relative homology on a sub-complex. The flag of each cell is stored in an appropriate structure and getters are provided to access to this information.
The Hdvf_core class stores the associated reduction in sparse matrices: row-major for \(f\), and column-major for \(g\), \(h\) and \(\partial'\). Getters are provided to access this information. However, according to the chosen HDVF computation option (OPT_BND, OPT_F, OPT_G, OPT_FULL) the reduction can be computed only partially (and thus faster).
The class provides perfect HDVF constuction operations: compute_perfect_hdvf() and compute_rand_perfect_hdvf(), which build perfect HDVFs by pairing iteratively critical cells through the A() operation.
If the user wishes to build an HDVF using other criteria, several find_pair_A() functions are provided (searching for valid pairs of cells for Arespecting various constraints). The A operation can be applied to any pair returned by these functions.
Homology/cohomology generators are actually algebraic objects, namely chains. Methods get_homology_chain() and get_cohomology_chain() return the homology and cohomology generator chain associated to a given critical cell. VTK export functions output all the cells of such chains with non zero coefficients.
HDVF | CoefficientType | a model of the Ring concept providing the ring used to compute homology. |
| ComplexType | a model of the AbstractChainComplex concept, providing the type of abstract chain complex used. |
| ChainType | a model of the SparseChain concept (by default, OSM::Sparse_chain), providing the type of sparse chains used (should be coherent with SparseMatrixType). |
| SparseMatrixType | a model of the SparseMatrix concept (by default, OSM::Sparse_matrix), providing the type of sparse matrices used. |
Public Types | |
| typedef ChainType< CoefficientType, OSM::COLUMN > | CChain |
| Type of column-major chains. | |
| typedef ChainType< CoefficientType, OSM::ROW > | RChain |
| Type of row-major chains. | |
| typedef SparseMatrixType< CoefficientType, OSM::COLUMN > | CMatrix |
| Type of column-major sparse matrices. | |
| typedef SparseMatrixType< CoefficientType, OSM::ROW > | RMatrix |
| Type of row-major sparse matrices. | |
Public Member Functions | |
| Hdvf_core (const ComplexType &K, int hdvf_opt=OPT_FULL) | |
| Default constructor. | |
| Hdvf_core (const Hdvf_core &hdvf) | |
| ~Hdvf_core () | |
| virtual PairCell | find_pair_A (int q, bool &found) const |
| Finds a valid PairCell of dimension q / q+1 for A. | |
| virtual PairCell | find_pair_A (int q, bool &found, size_t gamma) const |
Finds a valid PairCell for A containing gamma (a cell of dimension q) | |
| virtual std::vector< PairCell > | find_pairs_A (int q, bool &found) const |
| Finds all valid PairCell of dimension q / q+1 for A. | |
| virtual std::vector< PairCell > | find_pairs_A (int q, bool &found, size_t gamma) const |
Finds all valid PairCell for A containing gamma (a cell of dimension q) | |
| void | A (size_t gamma1, size_t gamma2, int q) |
| A operation: pairs critical cells. | |
| std::vector< PairCell > | compute_perfect_hdvf (bool verbose=false) |
| Computes a perfect HDVF. | |
| std::vector< PairCell > | compute_rand_perfect_hdvf (bool verbose=false) |
| Computes a random perfect HDVF. | |
| bool | is_perfect_hdvf () |
| Tests if a HDVF is perfect. | |
| virtual std::vector< std::vector< size_t > > | get_flag (FlagType flag) const |
Gets cells with a given flag in any dimension. | |
| virtual std::vector< size_t > | get_flag_dim (FlagType flag, int q) const |
Gets cells with a given flag in dimension q. | |
| FlagType | get_cell_flag (int q, size_t tau) const |
Gets the flag of the cell tau in dimension q. | |
| int | get_hdvf_opts () const |
| Gets HDVF computation option. | |
| const RMatrix & | get_f (int q) const |
| Gets the row-major matrix of \(f\) (from the reduction associated to the HDVF). | |
| const CMatrix & | get_g (int q) const |
| Gets the column-major matrix of \(g\) (from the reduction associated to the HDVF). | |
| const CMatrix & | get_h (int q) const |
| Gets the column-major matrix of \(h\) (from the reduction associated to the HDVF). | |
| const CMatrix & | get_dd (int q) const |
| Gets the column-major matrix of \(\partial'\), reduced boundary operator (from the reduction associated to the HDVF). | |
| std::ostream & | print_matrices (std::ostream &out=std::cout) const |
| Prints the matrices of the reduction. | |
| std::ostream & | print_reduction (std::ostream &out=std::cout) const |
| Prints the homology and cohomology reduction information. | |
| virtual std::vector< std::vector< int > > | get_psc_labels () const |
| Exports primary/secondary/critical labels (in particular for vtk export) | |
| virtual CChain | get_homology_chain (size_t cell, int q) const |
Gets homology generators associated to cell (critical cell) of dimension q (used by vtk export). | |
| virtual CChain | get_cohomology_chain (size_t cell, int dim) const |
Gets cohomology generators associated to cell (critical cell) of dimension q (used by vtk export). | |
| std::ostream & | save_hdvf_reduction (std::ostream &out) |
| Saves a HDVF together with the associated reduction (f, g, h, d matrices) | |
| std::istream & | load_hdvf_reduction (std::istream &out) |
| Loads a HDVF together with the associated reduction (f, g, h, d matrices) | |
Protected Member Functions | |
| template<int ChainTypeFlag> | |
| ChainType< CoefficientType, ChainTypeFlag > | projection (const ChainType< CoefficientType, ChainTypeFlag > &chain, FlagType flag, int q) const |
| void | progress_bar (size_t i, size_t n) |
Protected Attributes | |
| std::vector< std::vector< FlagType > > | _flag |
| std::vector< size_t > | _nb_P |
| std::vector< size_t > | _nb_S |
| std::vector< size_t > | _nb_C |
| std::vector< RMatrix > | _F_row |
| std::vector< CMatrix > | _G_col |
| std::vector< CMatrix > | _H_col |
| std::vector< CMatrix > | _DD_col |
| const ComplexType & | _K |
| int | _hdvf_opt |
| CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::Hdvf_core | ( | const ComplexType & | K, |
| int | hdvf_opt = OPT_FULL |
||
| ) |
Default constructor.
Builds an "empty" HDVF_core associated to K (with all cells critical). By default, the HDVF option is set to OPT_FULL (full reduction computed).
| [in] | K | A chain complex (a model of AbstractChainComplex) |
| [in] | hdvf_opt | Option for HDVF computation (OPT_BND, OPT_F, OPT_G or OPT_FULL) |
| void CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::A | ( | size_t | gamma1, |
| size_t | gamma2, | ||
| int | q | ||
| ) |
A operation: pairs critical cells.
A pair of critical cells \((\gamma_1, \gamma_2)\) of respective dimension q and q+1 is valid for A if \(\langle \partial_{q+1}(\gamma_2), \gamma_1 \rangle\) is invertible. After the A() operation, \(\gamma_1\) becomes PRIMARY, \(\gamma_2\) becomes SECONDARY. The A method updates the reduction accordingly (in time \(\mathcal O(n^2)\)).
| [in] | gamma1 | First cell of the pair (dimension q) |
| [in] | gamma2 | Second cell of the pair (dimension q+1) |
| [in] | q | Dimension of the pair |
| std::vector< PairCell > CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::compute_perfect_hdvf | ( | bool | verbose = false | ) |
Computes a perfect HDVF.
As long as valid pairs for A exist, the function selects the first available pair (returned by find_pair_A()) and applies the corresponding A() operation. If the Ring of coefficients is a field, this operation always produces a perfect HDVF (ie. the reduced boundary is null and the reduction provides homology and cohomology information). Otherwise the operation produces a maximal HDVF with a residual boundary matrix over critical cells.
If the HDVF is initially not trivial (some cells have already been paired), the function completes it into a perfect HDVF.
| [in] | verbose | If this parameter is true, all intermediate reductions are printed out. |
PairCell paired with A. | std::vector< PairCell > CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::compute_rand_perfect_hdvf | ( | bool | verbose = false | ) |
Computes a random perfect HDVF.
As long as valid pairs for A exist, the function selects a random pair (among pairs returned by find_pairs_A()) and applies the corresponding A() operation. If the Ring of coefficients is a field, this operation always produces a perfect HDVF (that is the reduced boundary is null and the reduction provides homology and cohomology information).
If the HDVF is initially not trivial (some cells have already been paired), the function randomly completes it into a perfect HDVF.
compute_perfect_hdvf() (finding out all possible valid pairs requires additional time).| [in] | verbose | If this parameter is true, all intermediate reductions are printed out. |
|
virtual |
Finds a valid PairCell of dimension q / q+1 for A.
The function searches a pair of critical cells \((\gamma_1, \gamma2)\) of dimension q / q+1, valid for A (ie. such that \(\langle \partial_{q+1}(\gamma_2), \gamma_1 \rangle\) invertible). It returns the first valid pair found by iterators.
| [in] | q | Lower dimension of the pair. |
| [in] | found | Reference to a Boolean variable. The method sets found to true if a valid pair is found, false otherwise. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
|
virtual |
Finds a valid PairCell for A containing gamma (a cell of dimension q)
The function searches a cell \(\gamma'\) such that one of the following conditions holds:
| [in] | q | Dimension of the cell gamma. |
| [in] | found | Reference to a Boolean variable. The method sets found to true if a valid pair is found, false otherwise. |
| [in] | gamma | Index of a cell to pair. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
|
virtual |
Finds all valid PairCell of dimension q / q+1 for A.
The function searches all pairs of critical cells \((\gamma_1, \gamma2)\) of dimension q / q+1, valid for A (ie. such that \(\langle \partial_{q+1}(\gamma_2), \gamma_1 \rangle\) invertible). It returns a vector of such pairs.
| [in] | q | Lower dimension of the pair. |
| [in] | found | Reference to a Boolean variable. The method sets found to true if a valid pair is found, false otherwise. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
|
virtual |
Finds all valid PairCell for A containing gamma (a cell of dimension q)
The function searches all CRITICAL cells \(\gamma'\) such that one of the following conditions holds:
| [in] | q | Dimension of the cell gamma. |
| [in] | found | Reference to a Boolean variable. The method sets found to true if a valid pair is found, false otherwise. |
| [in] | gamma | Index of a cell to pair. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
| FlagType CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::get_cell_flag | ( | int | q, |
| size_t | tau | ||
| ) | const |
Gets the flag of the cell tau in dimension q.
| [in] | tau | Index of the cell. |
| [in] | q | Dimension of the cell. |
|
virtual |
Gets cohomology generators associated to cell (critical cell) of dimension q (used by vtk export).
The method exports the chain \(f^\star(\sigma)\) for \(\sigma\) the cell of index cell and dimension q.
| [in] | cell | Index of the (critical) cell. |
| [in] | dim | Dimension of the (critical) cell. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >, and CGAL::HDVF::Hdvf_persistence< CoefficientType, ComplexType, DegType, FiltrationType >.
|
virtual |
Gets cells with a given flag in any dimension.
The function returns a vector containing, for each dimension, the vector of cells with a given flag.
| [in] | flag | Flag to select. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
|
virtual |
Gets cells with a given flag in dimension q.
The function returns the vector of cells of dimension q with a given flag.
| [in] | flag | Flag to select. |
| [in] | q | Dimension visited. |
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >.
|
virtual |
Gets homology generators associated to cell (critical cell) of dimension q (used by vtk export).
The method exports the chain \(g(\sigma)\) for \(\sigma\) the cell of index cell and dimension q.
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >, and CGAL::HDVF::Hdvf_persistence< CoefficientType, ComplexType, DegType, FiltrationType >.
|
virtual |
Exports primary/secondary/critical labels (in particular for vtk export)
The method exports the labels of every cells in each dimension.
Reimplemented in CGAL::HDVF::Hdvf_duality< CoefficientType, ComplexType >, and CGAL::HDVF::Hdvf_persistence< CoefficientType, ComplexType, DegType, FiltrationType >.
| bool CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::is_perfect_hdvf | ( | ) |
Tests if a HDVF is perfect.
The function returns true is the reduced boundary matrix is null and false otherwise.
| std::istream & CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::load_hdvf_reduction | ( | std::istream & | out | ) |
| std::ostream & CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::print_matrices | ( | std::ostream & | out = std::cout | ) | const |
Prints the matrices of the reduction.
Prints the matrices of the reduction (that is \(f\), \(g\), \(h\), \(\partial'\) the reduced boundary). By default, outputs the complex to std::cout.
| std::ostream & CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::print_reduction | ( | std::ostream & | out = std::cout | ) | const |
Prints the homology and cohomology reduction information.
Prints \(f^*\), \(g\) \(\partial'\) the reduced boundary over each critical cell.
By default, outputs the complex to std::cout.
| std::ostream & CGAL::HDVF::Hdvf_core< CoefficientType, ComplexType, ChainType, SparseMatrixType >::save_hdvf_reduction | ( | std::ostream & | out | ) |