CASM  1.1.0
A Clusters Approach to Statistical Mechanics
HermiteCounter.hh
Go to the documentation of this file.
1 #ifndef HERMITECOUNTER_HH
2 #define HERMITECOUNTER_HH
3 
5 #include "casm/external/Eigen/Dense"
7 
8 namespace CASM {
9 namespace xtal {
41  public:
42  typedef Eigen::VectorXi::Scalar value_type;
43  typedef CASM::Index Index;
44 
46  // HermiteCounter() {};
47 
50  HermiteCounter(int init_determinant, int init_dim);
51 
52  // You probably will never need these. They're just here for testing more than
53  // anything. Either way, they're safe to call.
54  Index position() const;
55  Eigen::VectorXi diagonal() const;
56  // value_type low() const;
57  // value_type high() const;
58 
60  Eigen::MatrixXi current() const;
61 
63  value_type determinant() const;
64 
66  Index dim() const;
67 
69  void reset_current();
70 
73  void next_determinant();
74 
77  void jump_to_determinant(value_type new_det);
78 
81 
83  Eigen::MatrixXi operator()() const;
84 
85  private:
89 
91  Eigen::VectorXi m_diagonal;
92 
96 
100 };
101 
102 namespace HermiteCounter_impl {
105 HermiteCounter::Index _spill_factor(Eigen::VectorXi &diag,
106  HermiteCounter::Index position,
108 
111 HermiteCounter::Index next_spill_position(Eigen::VectorXi &diag,
112  HermiteCounter::Index position);
113 
117 
120 EigenVectorXiCounter _upper_tri_counter(const Eigen::VectorXi &current_diag);
121 
124 Eigen::MatrixXi _zip_matrix(const Eigen::VectorXi &current_diag,
125  const Eigen::VectorXi &current_upper_tri);
126 
129 Eigen::MatrixXi _expand_dims_old(const Eigen::MatrixXi &hermit_mat,
130  const Eigen::VectorXi &active_dims);
131 
134 Eigen::MatrixXi _expand_dims(const Eigen::MatrixXi &H,
135  const Eigen::MatrixXi &G);
136 
139 Eigen::VectorXi _canonical_unroll(const Eigen::MatrixXi &hermit_mat);
140 
143 bool _canonical_compare(const Eigen::MatrixXi &H0, const Eigen::MatrixXi &H1);
144 } // namespace HermiteCounter_impl
145 
146 } // namespace xtal
147 } // namespace CASM
148 
149 #endif
A Counter allows looping over many incrementing variables in one loop.
Definition: Counter.hh:109
Index m_pos
Keeps track of the current diagonal element that needs to be factored.
Index _increment_diagonal()
Go to the next values of diagonal elements that keep the same determinant.
HermiteCounter(int init_determinant, int init_dim)
constructor to satisfy iterator requirements. Do not recommend.
HermiteCounter & operator++()
Jump to the next available HNF matrix.
void next_determinant()
Skip the remaining iterations and start at the next determinant value.
void jump_to_determinant(value_type new_det)
Reset the diagonal to the specified determinant and set the other elements to zero.
Eigen::MatrixXi operator()() const
Get the current matrix the counter is on.
EigenVectorXiCounter m_upper_tri
unrolled vector of the upper triangle (does not include diagonal elements)
Eigen::VectorXi::Scalar value_type
Index dim() const
Get the dimensions of *this.
Eigen::VectorXi m_diagonal
Vector holding diagonal element values.
Eigen::VectorXi diagonal() const
value_type determinant() const
Get the current determinant.
Eigen::MatrixXi current() const
Get the current matrix the counter is on.
void reset_current()
reset the counter to the first iteration of the current determinant
EigenVectorXiCounter _upper_tri_counter(const Eigen::VectorXi &current_diag)
Create a counter for the elements above the diagonal based on the current diagonal value.
Eigen::MatrixXi _expand_dims(const Eigen::MatrixXi &H, const Eigen::MatrixXi &G)
Expand a n x n Hermite normal matrix (H) into a m x m one through a m x m generating matrix (G) (e....
HermiteCounter::Index next_spill_position(Eigen::VectorXi &diag, HermiteCounter::Index position)
Spill the next factor of the specified element with its neighbor, and return new position.
bool _canonical_compare(const Eigen::MatrixXi &H0, const Eigen::MatrixXi &H1)
Compare two integer matrices and see which one is lexicographically greatest. Returns true if H0<H1.
Eigen::MatrixXi _zip_matrix(const Eigen::VectorXi &current_diag, const Eigen::VectorXi &current_upper_tri)
Assemble a matrix diagonal and unrolled upper triangle values into a matrix.
HermiteCounter::Index _spill_factor(Eigen::VectorXi &diag, HermiteCounter::Index position, HermiteCounter::value_type attempt)
Find the next factor of the specified position and share with next element. Use attempt as starting p...
Eigen::MatrixXi _expand_dims_old(const Eigen::MatrixXi &hermit_mat, const Eigen::VectorXi &active_dims)
Expand a n x n Hermite normal matrix into a m x m one (e.g. for 2D supercells)
Eigen::VectorXi _canonical_unroll(const Eigen::MatrixXi &hermit_mat)
Unroll a Hermit normal form square matrix into a vector such that it's canonical form is easy to comp...
HermiteCounter::Index upper_size(HermiteCounter::Index init_dim)
Determine the number of elements in the upper triangular matrix (excluding diagonal)
Main CASM namespace.
Definition: APICommand.hh:8
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39