StatMech
Loading...
Searching...
No Matches
ParityTransSector< RealType > Class Template Reference

Translation invariant sector of a many-body Hilbert space. More...

#include <ParityTransSector.hpp>

Inheritance diagram for ParityTransSector< RealType >:
Inheritance graph
Collaboration diagram for ParityTransSector< RealType >:
Collaboration graph

Public Types

using TotalSpace = ManyBodySpinSpace
 
- Public Types inherited from SubSpace< TotalSpace_, ScalarType_ >
using TotalSpace = TotalSpace_
 
using Scalar = ScalarType_
 
using Real = typename Eigen::NumTraits< ScalarType_ >::Real
 

Public Member Functions

__host__ __device__ ParityTransSector (int parity, ManyBodySpinSpace const &mbHSpace)
 Construct a ParityTransSector object from a ManyBodySpinSpace object.
 
__host__ __device__ ParityTransSector (int parity=0, int systemSize=0, int dimLoc=0)
 Default constructor.
 
__host__ __device__ constexpr int momentum () const
 
__host__ __device__ constexpr int parity () const
 
__host__ __device__ int rep (int j, int trans=0) const
 
__host__ __device__ int period (int j) const
 
- Public Member Functions inherited from SubSpace< TotalSpace_, ScalarType_ >
__host__ __device__ SubSpace (TotalSpace const &totalSpace)
 Construct a SubSpace object by copying totalSpace to its member variable (m_totalSpace)
 
__host__ __device__ SubSpace (TotalSpace &&totalSpace)
 Construct a SubSpace object by moving totalSpace to its member variable (m_totalSpace)
 
__host__ __device__ SubSpace ()=default
 Default constructor.
 
__host__ __device__ SubSpace (SubSpace const &other)=default
 Copy constructor.
 
template<typename Scalar_ >
__host__ __device__ SubSpace (SubSpace< TotalSpace_, Scalar_ > const &other)
 
__host__ __device__ SubSpace (SubSpace &&other)=default
 Move constructor.
 
__host__ __device__ ~SubSpace ()=default
 Destructor.
 
__host__ __device__ SubSpaceoperator= (SubSpace const &other)
 Copy assignment operator.
 
__host__ __device__ SubSpaceoperator= (SubSpace &&other)
 Move assignment operator.
 
__host__ __device__ bool operator== (SubSpace const &other) const
 Equality operator.
 
__host__ __device__ Matrix_tbasis ()
 
__host__ __device__ Matrix_t const & basis () const
 
__host__ __device__ TotalSpace const & totalSpace () const
 
__host__ __device__ int dimTot () const
 
- Public Member Functions inherited from HilbertSpace< Derived >
__host__ __device__ int dim () const
 
__host__ __device__ bool operator== (HilbertSpace const &other) const
 

Private Attributes

int m_parity
 
int m_numParityPairs = 0
 

Additional Inherited Members

- Protected Types inherited from SubSpace< TotalSpace_, ScalarType_ >
using Matrix_t = SparseCompressed< Scalar >
 
using Matrix_t = Eigen::SparseMatrix< Scalar >
 
- Protected Attributes inherited from SubSpace< TotalSpace_, ScalarType_ >
TotalSpace m_totalSpace
 
Matrix_t m_basisStates
 

Detailed Description

template<typename RealType = Real_t>
class ParityTransSector< RealType >

Translation invariant sector of a many-body Hilbert space.

Template Parameters
Scalar_t
int
Matrix

Member Typedef Documentation

◆ TotalSpace

template<typename RealType = Real_t>
using ParityTransSector< RealType >::TotalSpace = ManyBodySpinSpace

Constructor & Destructor Documentation

◆ ParityTransSector() [1/2]

template<typename RealType >
__host__ __device__ ParityTransSector< RealType >::ParityTransSector ( int  parity,
ManyBodySpinSpace const &  mbHSpace 
)

Construct a ParityTransSector object from a ManyBodySpinSpace object.

Parameters
k
mbHSpace
65 m_parity((parity >= 0 ? +1 : -1)) {
66 debug_constructor_printf(1);
67 constexpr Complex_t<RealType> I = Complex_t<RealType>(0, 1);
68 if(this->totalSpace().sysSize() == 0) { return; }
69
71 debug_printf("%s:\n\t this->totalSpace().dim()=%d, this->totalSpace().transEqDim()=%d, "
72 "this->totalSpace().transPeriod().sum()=%d\n",
73 __PRETTY_FUNCTION__, this->totalSpace().dim(), this->totalSpace().transEqDim(),
74 this->totalSpace().transPeriod().sum());
75 assert(this->totalSpace().transPeriod().sum() == this->totalSpace().dim());
76
77 int L = this->totalSpace().sysSize();
78 int dim = 0;
80 Eigen::ArrayX<bool> calculated;
81 calculated = Eigen::ArrayX<bool>::Constant(this->totalSpace().transEqDim(), false);
82 for(int eqClass = 0; eqClass < this->totalSpace().transEqDim(); ++eqClass) {
83 if(calculated(eqClass)) continue;
84 calculated(eqClass) = true;
85 auto repState = this->totalSpace().transEqClassRep(eqClass);
86 auto revState = this->totalSpace().reverse(repState);
87 auto reversed = this->totalSpace().stateToTransEqClass(revState);
88 calculated(reversed) = true;
89 // std::cout << "reversed = " << reversed << std::endl;
90 // std::cout << this->totalSpace().ordinalToConfig(repState) << std::endl;
91 // std::cout << this->totalSpace().ordinalToConfig(revState) << "\n" << std::endl;
92 if(eqClass != reversed) { m_numParityPairs += 1; }
93 dim += 1;
94 }
95 if(m_parity == -1) dim = m_numParityPairs;
96 std::cout << "m_parity = " << m_parity << ", dim = " << dim
97 << ", m_numParityPairs = " << m_numParityPairs << std::endl;
98
99#ifdef __CUDA_ARCH__
100 this->basis().resize(this->totalSpace().dim(), dim, this->totalSpace().dim());
101
102 int nonZeroId = 0, classIdx = 0;
103 if(m_parity != -1) {
104 calculated = Eigen::ArrayX<bool>::Constant(this->totalSpace().transEqDim(), false);
105 for(int eqClass = 0; eqClass != this->totalSpace().transEqDim(); ++eqClass) {
106 if(calculated(eqClass)) continue;
107 calculated(eqClass) = true;
108 auto repState = this->totalSpace().transEqClassRep(eqClass);
109 auto revState = this->totalSpace().reverse(repState);
110 auto reversed = this->totalSpace().stateToTransEqClass(revState);
111 calculated(reversed) = true;
112 if(eqClass != reversed) continue;
113
114 this->basis().outerIndexPtr()[classIdx] = nonZeroId;
115 for(int trans = 0; trans != this->totalSpace().transPeriod(eqClass); ++trans) {
116 this->basis().innerIndexPtr()[nonZeroId]
117 = this->totalSpace().transEqClassRep(eqClass, trans);
118 this->basis().valuePtr()[nonZeroId]
119 = 1.0 / RealType(sqrt(RealType(this->totalSpace().transPeriod(eqClass))));
120 ++nonZeroId;
121 }
122 ++classIdx;
123 }
124 }
125
126 calculated = Eigen::ArrayX<bool>::Constant(this->totalSpace().transEqDim(), false);
127 for(int eqClass = 0; eqClass != this->totalSpace().transEqDim(); ++eqClass) {
128 if(calculated(eqClass)) continue;
129 calculated(eqClass) = true;
130 auto repState = this->totalSpace().transEqClassRep(eqClass);
131 auto revState = this->totalSpace().reverse(repState);
132 auto reversed = this->totalSpace().stateToTransEqClass(revState);
133 calculated(reversed) = true;
134 if(eqClass == reversed) continue;
135
136 this->basis().outerIndexPtr()[classIdx] = nonZeroId;
137 for(int trans = 0; trans != this->totalSpace().transPeriod(eqClass); ++trans) {
138 this->basis().innerIndexPtr()[nonZeroId]
139 = this->totalSpace().transEqClassRep(eqClass, trans);
140 this->basis().valuePtr()[nonZeroId]
141 = 1.0 / RealType(sqrt(RealType(2 * this->totalSpace().transPeriod(eqClass))));
142 ++nonZeroId;
143 }
144 for(int trans = 0; trans != this->totalSpace().transPeriod(eqClass); ++trans) {
145 this->basis().innerIndexPtr()[nonZeroId]
146 = this->totalSpace().transEqClassRep(reversed, trans);
147 this->basis().valuePtr()[nonZeroId]
148 = m_parity / RealType(sqrt(RealType(2 * this->totalSpace().transPeriod(reversed))));
149 ++nonZeroId;
150 }
151 ++classIdx;
152 }
153 assert(dim == classIdx && "dim == classIdx");
154 this->basis().outerIndexPtr()[dim] = nonZeroId;
155 this->basis().setNonZeros(this->basis().outerIndexPtr()[dim]);
156 for(int eqClass = 0; eqClass != dim; ++eqClass) {
157 thrust::sort_by_key(
158 thrust::seq, this->basis().innerIndexPtr() + this->basis().outerIndexPtr()[eqClass],
159 this->basis().innerIndexPtr() + this->basis().outerIndexPtr()[eqClass + 1],
160 this->basis().valuePtr() + this->basis().outerIndexPtr()[eqClass]);
161 }
162#else // #ifdef __CUDA_ARCH__
163 this->basis().resize(this->totalSpace().dim(), dim);
164 this->basis().reserve(Eigen::VectorXi::Constant(dim, 2 * L));
165
166 int classIdx = 0;
167 if(m_parity != -1) {
168 calculated = Eigen::ArrayX<bool>::Constant(this->totalSpace().transEqDim(), false);
169 for(int eqClass = 0; eqClass != this->totalSpace().transEqDim(); ++eqClass) {
170 if(calculated(eqClass)) continue;
171 calculated(eqClass) = true;
172 auto repState = this->totalSpace().transEqClassRep(eqClass);
173 auto revState = this->totalSpace().reverse(repState);
174 auto reversed = this->totalSpace().stateToTransEqClass(revState);
175 calculated(reversed) = true;
176 if(eqClass != reversed) continue;
177
178 for(int trans = 0; trans != this->totalSpace().transPeriod(eqClass); ++trans) {
179 auto state = this->totalSpace().transEqClassRep(eqClass, trans);
180 this->basis().insert(state, classIdx)
181 = 1.0 / RealType(sqrt(RealType(this->totalSpace().transPeriod(eqClass))));
182 }
183 ++classIdx;
184 }
185 }
186
187 calculated = Eigen::ArrayX<bool>::Constant(this->totalSpace().transEqDim(), false);
188 for(int eqClass = 0; eqClass != this->totalSpace().transEqDim(); ++eqClass) {
189 if(calculated(eqClass)) continue;
190 calculated(eqClass) = true;
191 auto repState = this->totalSpace().transEqClassRep(eqClass);
192 auto revState = this->totalSpace().reverse(repState);
193 auto reversed = this->totalSpace().stateToTransEqClass(revState);
194 calculated(reversed) = true;
195 if(eqClass == reversed) continue;
196
197 for(int trans = 0; trans != this->totalSpace().transPeriod(eqClass); ++trans) {
198 auto state = this->totalSpace().transEqClassRep(eqClass, trans);
199 this->basis().insert(state, classIdx)
200 = 1.0 / RealType(sqrt(RealType(2 * this->totalSpace().transPeriod(eqClass))));
201 }
202 for(int trans = 0; trans != this->totalSpace().transPeriod(reversed); ++trans) {
203 auto state = this->totalSpace().transEqClassRep(reversed, trans);
204 this->basis().insert(state, classIdx)
205 = m_parity / RealType(sqrt(RealType(2 * this->totalSpace().transPeriod(reversed))));
206 }
207 ++classIdx;
208 }
209
210 this->basis().makeCompressed();
211#endif // #ifdef __CUDA_ARCH__
212}
Definition mytypes.hpp:147
__host__ __device__ int dim() const
Definition HilbertSpace.hpp:34
__host__ __device__ Vector_t const & stateToTransEqClass() const
Definition HilbertSpace.hpp:320
__host__ __device__ Vector_t const & transPeriod() const
Definition HilbertSpace.hpp:317
__host__ __device__ int transEqDim() const
Definition HilbertSpace.hpp:309
__host__ __device__ void computeTransEqClass() const
Definition HilbertSpace.hpp:357
__host__ __device__ int reverse(int state) const
Reverse the input state.
Definition HilbertSpace.hpp:347
__host__ __device__ int sysSize() const
Definition HilbertSpace.hpp:265
__host__ __device__ Vector_t const & transEqClassRep() const
Definition HilbertSpace.hpp:311
int m_parity
Definition ParityTransSector.hpp:21
int m_numParityPairs
Definition ParityTransSector.hpp:22
__host__ __device__ constexpr int parity() const
Definition ParityTransSector.hpp:47
__host__ __device__ int * innerIndexPtr() const
Definition MatrixUtils.cuh:421
__host__ __device__ void resize(int rows, int cols, int reserved=1)
Definition MatrixUtils.cuh:353
__host__ __device__ int & setNonZeros(int input)
Definition MatrixUtils.cuh:417
__host__ __device__ int * outerIndexPtr() const
Definition MatrixUtils.cuh:420
__host__ __device__ Scalar_t * valuePtr() const
Definition MatrixUtils.cuh:422
Definition HilbertSpace.hpp:568
__host__ __device__ Matrix_t & basis()
Definition HilbertSpace.hpp:668
__host__ __device__ TotalSpace const & totalSpace() const
Definition HilbertSpace.hpp:671

◆ ParityTransSector() [2/2]

template<typename RealType = Real_t>
__host__ __device__ ParityTransSector< RealType >::ParityTransSector ( int  parity = 0,
int  systemSize = 0,
int  dimLoc = 0 
)
inline

Default constructor.

Parameters
k
systemSize
dimLoc
41 : ParityTransSector(parity, ManyBodySpinSpace(systemSize, dimLoc)) {
42 debug_constructor_printf((Default));
43 }
Definition HilbertSpace.hpp:423
Translation invariant sector of a many-body Hilbert space.
Definition ParityTransSector.hpp:19

Member Function Documentation

◆ momentum()

template<typename RealType = Real_t>
__host__ __device__ constexpr int ParityTransSector< RealType >::momentum ( ) const
inlineconstexpr
45{ return 0; };

◆ parity()

template<typename RealType = Real_t>
__host__ __device__ constexpr int ParityTransSector< RealType >::parity ( ) const
inlineconstexpr
47{ return m_parity; };

◆ period()

template<typename RealType = Real_t>
__host__ __device__ int ParityTransSector< RealType >::period ( int  j) const
inline
55 {
56 assert(j < this->dim() && "j < this->dim()");
57 return this->basis().outerIndexPtr()[j + 1] - this->basis().outerIndexPtr()[j];
58 };

◆ rep()

template<typename RealType = Real_t>
__host__ __device__ int ParityTransSector< RealType >::rep ( int  j,
int  trans = 0 
) const
inline
49 {
50 int innerId = this->basis().outerIndexPtr()[j] + (trans % this->period(j));
51 assert(innerId < this->basis().nonZeros() && "innerId < this->basis().nonZeros()");
52 return this->basis().innerIndexPtr()[innerId];
53 };
__host__ __device__ int period(int j) const
Definition ParityTransSector.hpp:55

Member Data Documentation

◆ m_numParityPairs

template<typename RealType = Real_t>
int ParityTransSector< RealType >::m_numParityPairs = 0
private

◆ m_parity

template<typename RealType = Real_t>
int ParityTransSector< RealType >::m_parity
private

The documentation for this class was generated from the following file: