Skip to content

Electron

A block in an Electron (fermions with \(\uparrow, \downarrow\) spin) Hilbert space.

Sources
electron.hpp
electron.cpp
electron.jl


Constructors

Electron(int64_t nsites, std::string backend = "auto");
Electron(int64_t nsites, int64_t nup, int64_t ndn, std::string backend = "auto");
Electron(int64_t nsites, Representation irrep, std::string backend = "auto");
Electron(int64_t nsites, int64_t nup, int64_t ndn, Representation irrep, std::string backend = "auto");
Electron(nsites::Int64, backend::String="auto")
Electron(nsites::Int64, nup::Int64, ndn::Int64, backend::String="auto")
Electron(nsites::Int64, irrep::Representation, backend::String="auto")
Electron(nsites::Int64, nup::Int64, ndn::Int64, irrep::Representation, backend::String="auto")
Name Description
nsites number of sites (integer)
nup number of "up" electrons (integer)
ndn number of "dn" electrons (integer)
irrep Irreducible Representation of the symmetry group
backend backend used for coding the basis states auto

The parameter backend chooses how the block is coded internally. By using the default parameter auto the backend is chosen automatically. Alternatives are 32bit, 64bit


Iteration

An Electron block can be iterated over, where at each iteration a ProductState representing the corresponding basis state is returned.

auto block = Electron(4, 2, 2);
for (auto pstate : block) {
  Log("{} {}", to_string(pstate), index(block, pstate));
}
block = Electron(4, 2, 2)
for pstate in block
    @show pstate, index(block, pstate) 
end

Methods

index

Returns the index of a given ProductState in the basis of the tJ block.

int64_t index(tJ const &block, ProductState const &pstate);
index(block::tJ, pstate::ProductState)::Int64

1-indexing

In the C++ version, the index count starts from "0" whereas in Julia the index count starts from "1".


nsites

Returns the number of sites of the block.

int64_t nsites(tJ const &block);
nsites(block::tJ)::Int64

size

Returns the size of the block, i.e. its dimension.

int64_t size(tJ const &block) const;
size(block::tJ)::Int64

dim

Returns the dimension of the block, same as "size" for non-distributed blocks.

int64_t dim(tJ const &block) const;
dim(block::tJ)::Int64

isreal

Returns whether the block can be used with real arithmetic. Complex arithmetic is needed when a Representation is genuinely complex.

bool isreal(tJ const &block);
isreal(block::tJ)::Bool

Usage Example

int N = 4;
int nup = 2;
int ndn = 1;

// without number conservation
auto block = Electron(N);
XDIAG_SHOW(block);

// with number conservation
auto block_np = Electron(N, nup, ndn);
XDIAG_SHOW(block_np);

// with symmetries, without number conservation
auto p1 = Permutation({0, 1, 2, 3});
auto p2 = Permutation({1, 2, 3, 0});
auto p3 = Permutation({2, 3, 0, 1});
auto p4 = Permutation({3, 0, 1, 2});
auto group = PermutationGroup({p1, p2, p3, p4});
auto irrep = Representation(group, arma::vec{1, -1, 1, -1});
auto block_sym = Electron(N, irrep);
XDIAG_SHOW(block_sym);

// with symmetries and number conservation
auto block_sym_np = Electron(N, nup, ndn, irrep);
XDIAG_SHOW(block_sym_np);
XDIAG_SHOW(block_sym_np.nsites());
XDIAG_SHOW(block_sym_np.size());

// Iteration
for (auto pstate : block_sym_np) {
  Log("{} {}", to_string(pstate), index(block_sym_np, pstate));
}
N = 4
nup = 2
ndn = 1

# without number conservation
block = Electron(N)
@show block

# with number conservation
block_np = Electron(N, nup, ndn)
@show block_np

# with symmetries, without number conservation
p = Permutation([2, 3, 4, 1])
group = PermutationGroup([p^0, p^1, p^2, p^3])
rep = Representation(group, [1.0, -1.0, 1.0, -1.0])
block_sym = Electron(N, rep)
@show block_sym

# with symmetries and number conservation
block_sym_np = Electron(N, nup, ndn, rep)
@show block_sym_np
@show nsites(block_sym_np)
@show size(block_sym_np)

# Iteration
for pstate in block_sym_np
    @show pstate, index(block_sym_np, pstate)
end