Complex#
Complex#
-
class Complex#
Primary class for computing persistent Laplacians
- Template Parameters:
eigs_Algorithm – Algorithm wrapper class to use for computing eigenvalues of the Laplacian in the spectra family of functions. Default is a wrapper for Eigen::SelfAdjointEigenSolver.
up_Algorithm – Algorithm wrapper class to use for computing the up-Laplacian. Default is the Schur complement algorithm presented in Memoli, Wan, and Wang 2020. //TODO: \see up_algorithms and eigs_algorithms
Subclassed by petls::Alpha, petls::Rips, petls::dFlag
Public Functions
-
Complex()#
Default constructor with no boundary maps or simplices.
-
Complex(std::vector<Eigen::SparseMatrix<storage>> boundaries, std::vector<std::vector<filtration_type>> filtrations)#
Primary constructor. Important Assumptions: 1) Boundary matrix has real coefficients stored as integers (but not mod 2!) 2) Boundary matrix dimensions agree with filtrations sizes 3) Length(filtrations) = Length(boundaries) + 1
- Parameters:
boundaries – a vector of Eigen::SparseMatrix of type int. Boundaries must be sorted in order of dimension.
filtrations – a vector of vector of filtrations. filtrations[dim] is a list of all filtrations of simplices in dimension dim. filtrations must be sorted in order of dimension, each filtrations[i] must be sorted in order of filtration
-
void set_boundaries_filtrations(std::vector<Eigen::SparseMatrix<storage>> boundaries, std::vector<std::vector<filtration_type>> filtrations)#
Set the boundaries and filtrations of a complex, particularly if the default constructor was called.
- Parameters:
boundaries – a vector of Eigen::SparseMatrix of type int. Boundaries must be sorted in order of dimension.
filtrations – a vector of vector of filtrations. filtrations[dim] is a list of all filtrations of simplices in dimension dim. filtrations must be sorted in order of dimension, each filtrations[i] must be sorted in order of filtration. Important Assumptions: 1) Boundary matrix has real coefficients stored as integers (but not mod 2!), 2) Boundary matrix dimensions agree with filtrations sizes, 3) Length(filtrations) = Length(boundaries) + 1.
-
inline void set_verbose(bool verbose)#
Set verbose
- Parameters:
verbose – New setting.
-
inline void set_flipped(bool use_flipped)#
Set flipped
- Parameters:
use_flipped – New setting
-
void get_L(int dim, filtration_type a, filtration_type b, DenseMatrix_PL &L)#
Get the Persistent Laplacian Matrix (by reference).
- Parameters:
dim – dimension
a – start filtration value
b – end filtration value (must be >= a)
L – [out] the matrix where the persistent Laplacian will be stored (by reference).
-
DenseMatrix_PL get_L(int dim, filtration_type a, filtration_type b)#
Get the Persistent Laplacian Matrix (by value).
Warning: this does a potentially expensive copy. The primary reason this function exists is to provide reasonable python binding access to the matrix L itself. If you do not need the persistent Laplacian matrix directly in python (e.g. you want its eigenvalues) it will be more efficient to compute by reference via get_L(int dim, filtration_type a, filtration_type b, DenseMatrix_PL &L).
- Parameters:
dim – dimension
a – start filtration value
b – end filtration value (must be >= a)
- Returns:
L the persistent Laplacian matrix.
-
void get_L_top_dim_flipped(filtration_type a, SparseMatrix_PL &L)#
Get a matrix with the same nonzero eigenvalues as the top-dimensional Persistent Laplacian Matrix (by reference).
- Parameters:
a – start filtration value
L – [out] the matrix (by reference) that has the same eigenvalues as the top-dimensional Persistent Laplacian
-
SparseMatrix_PL get_L_top_dim_flipped(filtration_type a)#
Get a matrix with the same nonzero eigenvalues as the top-dimensional Persistent Laplacian Matrix (by value).
Warning: this does a potentially expensive copy. The primary reason this function exists is to provide reasonable python binding access to the matrix itself. If you do not need this matrix directly in python (e.g. you want its eigenvalues) it will be more efficient to compute by reference via get_L_top_dim_flipped(filtration_type a, SparseMatrix_PL &L).
- Parameters:
a – start filtration value
- Returns:
A matrix with the same nonzero eigenvalues as the top-dimensional persistent Laplacian (but not the persistent Laplacian itself).
-
void get_up(int dim, filtration_type a, filtration_type b, DenseMatrix_PL &L_up)#
Get the up persistent Laplacian (by reference). The algorithm used is determined by the template parameter. See up_algorithms.hpp.
- Parameters:
dim – dimension
a – start filtration value
b – end filtration value
L_up – [out] up persistent Laplacian (by reference)
-
DenseMatrix_PL get_up(int dim, filtration_type a, filtration_type b)#
Get the up persistent Laplacian (by value). The algorithm used is determined by the template parameter. See up_algorithms.hpp.
- Parameters:
dim – dimension
a – start filtration value
b – end filtration value
- Returns:
up persistent Laplacian
-
void get_down(int dim, filtration_type a, Eigen::SparseMatrix<storage> &L_down)#
Get the down persistent Laplacian (by reference).
- Parameters:
dim – dimension
a – start filtration value
L_down – [out] the matrix (by reference) of the down presistent Laplacian
-
DenseMatrix_PL get_down(int dim, filtration_type a)#
Get the down persistent Laplacian (by value).
Warning: this does a potentially expensive copy. The primary reason this function exists is to provide reasonable python binding access to the matrix itself. If you do not need this matrix directly in python (e.g. you want its eigenvalues) it will be more efficient to compute by reference via get_down(int dim, filtration_type a, Eigen::SparseMatrix<storage> &L_down).
- Parameters:
dim – dimension
a – start filtration value \return The matrix (by value) of the down presistent Laplacian
-
std::vector<spectra_type> nonzero_spectra(int dim, filtration_type a, filtration_type b, SparseMatrix_PL PH_basis, bool use_dummy_harmonic_basis)#
Compute the nonzero eigenvalues of a persistent Laplacian using Schur restriction with the null space, either given by Persistent Homology representatives or by computing the null space.
- Parameters:
dim – dimension
a – start filtration level
b – end filtration level
(optional) – PH_basis Basis for the null space of the Laplacian, possibly obtained through persistent homology
use_dummy_harmonic_basis – Compute the null space of the Laplacian here, then do the same projected as if we had been given the null space
- Returns:
sorted vector of real, positive eigenvalues
-
std::vector<spectra_type> spectra(int dim, filtration_type a, filtration_type b)#
Get the Persistent Laplacian’s eigenvalues at a given dimension and filtration.
- Parameters:
dim – dimension
a – start filtration level
b – end filtration level
- Returns:
sorted vector of real, nonnegative eigenvalues
-
std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>>> spectra()#
Get all eigenvalues for all combinations of dimension and successive filtration values: a=filtrations[i] and b=filtrations[i+1]. Note: the caller does not know what spectra to expect from this.
- Returns:
vector of tuples (dim, a, b, eigenvalues) where “eigenvalues” is a sorted vector of real, nonnegative eigenvalues.
-
std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>>> spectra(std::vector<std::tuple<int, filtration_type, filtration_type>> spectra_request_list)#
This function essentially just calls spectra(dim, a, b) in a loop.
- Parameters:
spectra_quest_list – vector of tuples (dim, a, b) to compute the eigenvalues of L_{dim}^{a,b}.
- Returns:
vector of tuples (dim, a, b, eigenvalues), where eigenvalues it istelf a vector of real, nonnegative eigenvalues.
-
std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>>> spectra_allpairs()#
Get all eigenvalues for all combinations of dimension and filtration values. Note: the caller does not know what spectra to expect from this.
- Returns:
vector of tuples (dim, a, b, eigenvalues) where “eigenvalues” is a sorted vector of real, nonnegative eigenvalues.
-
std::pair<std::vector<spectra_type>, DenseMatrix_PL> eigenpairs(int dim, filtration_type a, filtration_type b)#
Get the Persistent Laplacian’s eigenvalues and eigenvectors at a given dimension and filtration.
- Parameters:
dim – dimension
a – start filtration level
b – end filtration level
- Returns:
sorted pair of: vector of real, nonnegative eigenvalues and Eigen matrix where column i is the eigenvector for eigenvalue i
-
std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>, DenseMatrix_PL>> eigenpairs()#
Get all eigenvalues for all combinations of dimension and successive filtration values: a=filtrations[i] and b=filtrations[i+1]. Note: the caller does not know what spectra to expect from this.
- Returns:
vector of tuples (dim, a, b, eigenvalues, eigenvectors) where “eigenvalues” is a sorted vector of real, nonnegative eigenvalues and “eigenvectors” is an Eigen::MatrixXf where column i is the eigenvector for eigenvalue i
-
std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>, DenseMatrix_PL>> eigenpairs(std::vector<std::tuple<int, filtration_type, filtration_type>> spectra_request_list)#
This function essentially just calls eigenpairs(dim, a, b) in a loop.
- Parameters:
spectra_quest_list – vector of tuples (dim, a, b) to compute the eigenvalues of L_{dim}^{a,b}.
- Returns:
vector of tuples (dim, a, b, eigenvalues, eigenvalues), where eigenvalues it istelf a vector of real, nonnegative eigenvalues and “eigenvectors” is an Eigen::MatrixXf where column i is the eigenvector for eigenvalue i.
-
std::pair<int, spectra_type> eigenvalues_summarize(std::vector<spectra_type> eigenvalues)#
Utility function to get the Betti number and least nonzero eigenvalue form a vector of eigenvalues.
- Parameters:
eigenvalues – eigenvalues
- Returns:
pair of (Betti number, least nonzero eigenvalue)
-
void store_L(int dim, filtration_type a, filtration_type b, std::string filename)#
Compute and store a persistent Laplacian matrix in a file in matrix market format.
- Parameters:
dim – dimension
a – start filtration value
b – end filtration value
filename – file to store matrix (typically .mtx extension)
-
void print_boundaries()#
Print all boundaries and corresponding filtrations.
-
void store_spectra(std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>>> spectra, std::string out_prefix)#
Write spectra to files.
for each dimension dim of the complex, files write “{out_prefix}_spectra_{dim}.txt”. Each line of the file is a space-separated list of eigenvalues. Lines may be empty. Note the filtration values are not reported.
- Parameters:
spectra – tuples (dim, a, b, eigenvalues)
out_prefix – Eigenvalues will be written to “{out_prefix}_spectra_{dim}.txt”
-
void store_spectra_summary(std::vector<std::tuple<int, filtration_type, filtration_type, std::vector<spectra_type>>> spectra, std::string out_prefix)#
Write spectra summary to files “{out_prefix}_spectra_summary.txt”
Each line is a space-separated list of filtrations, bettti numbers, and least nonzero eigenvalues: (filtration a) (filtration b) (betti 0) … (betti top_dim) (lambda 0) … (lambda top_dim)
- Parameters:
spectra – tuples (dim, a, b, eigenvalues)
out_prefix – Eigenvalues will be written to “{out_prefix}_spectra_summary.txt”
-
inline void time_to_csv(std::string filename)#
Write the profile to a csv file. This is done automatically when spectra() is called, but it can also be done here.BDCSVDEigen
- Parameters:
filename – name of the file to write to
-
std::vector<std::tuple<int, filtration_type, filtration_type>> filtration_list_to_spectra_request(std::vector<filtration_type> filtrations, std::vector<int> dims)#
Get tuples (dim, a, b) for all combinations of dimension and successive filtration values: a=filtrations[i], b=filtrations[i+1].
- Parameters:
filtrations – vector of filtration values
dims – vector of dimensions
- Returns:
tuples (dim, a, b) for all combinations of dimension and successive filtration values: a=filtrations[i], b=filtrations[i+1].
-
std::vector<std::tuple<int, filtration_type, filtration_type>> filtration_list_to_spectra_request_allpairs(std::vector<filtration_type> filtrations, std::vector<int> dims)#
Get tuples (dim, a, b) for all combinations of dimension and filtrations.
- Parameters:
filtrations – vector of filtration values
dims – vector of dimensions
- Returns:
tuples (dim, a, b) for all combinations of dimension and filtration values.
-
std::vector<filtration_type> get_all_filtrations()#
Get all unique filtration values in the complex.
Public Members
-
int top_dim#
Top dimension of the complex.
-
std::vector<FilteredBoundaryMatrix<storage>> filtered_boundaries#
Boundary matrix assuming real (or integer) coefficients.
-
bool use_flipped#
Compute the top-dimensional Laplacian’s eigenvalues in spectra function via the eigenvalues of the smaller of B_N B_N^T or B_N^T B_N and possible zero-padding.
-
Profile profile#
Profiler to track time usage of various steps.