Mackey  V3.3
A C++ library for computing RO(G) graded homology
Homology.hpp
Go to the documentation of this file.
1 #pragma once
3 #include "Utility/General.hpp"
4 #include "Smith.hpp"
5 #include "Abelian.hpp"
6 
9 
10 namespace mackey {
11 
13  //
18  template<typename rank_t, typename diff_t>
19  class Homology {
20 
21  typedef std::conditional_t<SFINAE::is_finite_cyclic<scalar_t<diff_t>>::value, scalar_t<diff_t>, int64_t> HScalar; //scalar with extended accuracy for Z coefficients (no change for Z/N coefficients)
22  typedef std::conditional_t<SFINAE::is_Dense<diff_t>::value, Eigen::Matrix<HScalar, -1, -1>, Eigen::SparseMatrix<HScalar, 0, storage_t<diff_t>>> diff_t_C; //diff_t column major with the scalar above
23  typedef std::conditional_t<SFINAE::is_Dense<diff_t>::value, Eigen::Matrix<HScalar, -1, -1, 1>, Eigen::SparseMatrix<HScalar, 1, storage_t<diff_t>>> diff_t_R; //diff_t row major with the scalar above
24  public:
25 
27  typedef diff_t_C Gens_t;
30 
33  bool isZero;
34 
36  Homology()=default;
37 
41  Homology(const Junction<rank_t, diff_t>& J, bool getQ=0);
42 
45 
48  rank_t basis(const gen_t&) const;
49 
51  gen_t boundary(const gen_t&) const;
52 
53  private:
54  diff_t_R Out_Qi, In_P_full, In_P_reduced;
55  std::vector<typename diff_t::StorageIndex> dontModOut;
56  diff_t_C In_Q;
57  row_vector_t<diff_t_C> diagonal;
58  typename diff_t::StorageIndex M;
59  diff_t_C getKernel(diff_t_C&);
60  void KernelModImage(diff_t_C&, diff_t_C&, bool);
61 
62  };
63 }
64 #include "impl/Homology.ipp"
Contains the class mackey::AbelianGroup.
Contains the classes mackey::Arrow, mackey::Chains and mackey::Junction.
Contains general functions and vector overloads independent of everything else.
Contains the class mackey::SmithNormalForm.
The Homology of a Junction.
Definition: Homology.hpp:19
AbelianGroup< rank_t > group
Encodes the homology groups as follows: group=[1,2,3] means homology Z+Z/2+Z/3. This works even for Z...
Definition: Homology.hpp:31
rank_t basis(const gen_t &) const
Given an element in homology, write it as a linear combination of the generators of the homology.
diff_t_C Gens_t
The type of our matrix of generators.
Definition: Homology.hpp:27
bool isZero
1 if the homology is trivial
Definition: Homology.hpp:33
Homology()=default
Default constructor.
gen_t boundary(const gen_t &) const
Given an x that is a boundary returns a y s.t. dy=x.
Homology(const Junction< rank_t, diff_t > &J, bool getQ=0)
Gens_t generators
Encodes the generators homology groups as follows: The i-th column corresponds to the generator for g...
Definition: Homology.hpp:32
col_vector_t< diff_t_C > gen_t
The dense type of our generators (a column in the generator matrix, always dense for convenience)
Definition: Homology.hpp:29
Everything in this library is under this namespace.
Definition: Box.hpp:9
Eigen::Matrix< scalar_t< T >, -1, 1 > col_vector_t
Dense column matrix.
Definition: Aliases.hpp:26
typename T::Scalar scalar_t
Scalar of matrix.
Definition: Aliases.hpp:14
Eigen::Matrix< scalar_t< T >, 1,-1 > row_vector_t
Dense row matrix.
Definition: Aliases.hpp:22