Skip to content

Electron

Representation of a block in a Electron (spinful fermion) Hilbert space.

Source electron.hpp

Constructors

Electron(n_sites::Integer)
Electron(n_sites::Integer, n_up::Integer, n_dn::Integer)
Electron(n_sites::Integer, group::PermutationGroup, irrep::Representation)
Electron(n_sites::Integer, n_up::Integer, n_dn::Integer, 
         group::PermutationGroup, irrep::Representation)
Electron(int64_t n_sites);
Electron(int64_t n_sites, int64_t n_up, int64_t n_dn);
Electron(int64_t n_sites, PermutationGroup permutation_group,
         Representation irrep);
Electron(int64_t n_sites, int64_t n_up, int64_t n_dn, 
         PermutationGroup group, Representation irrep);
Name Description
n_sites number of sites (integer)
n_up number of "up" electrons (integer)
n_dn number of "dn" electrons (integer)
group PermutationGroup defining the permutation symmetries
irrep Irreducible Representation of the symmetry group

Iteration

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

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

Methods

index

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

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

1-indexing

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

n_sites

Returns the number of sites of the block.

n_sites(block::Electron)
int64_t n_sites() const;

n_up

Returns the number of "up" electrons.

n_up(block::Electron)
int64_t n_up() const;

n_dn

Returns the number of "down" electrons.

n_dn(block::Electron)
int64_t n_dn() const;

permutation_group

Returns the PermutationGroup of the block, if defined.

permutation_group(block::Electron)
PermutationGroup permutation_group() const;

irrep

Returns the Representation of the block, if defined.

irrep(block::Electron)
Representation irrep() const;

size

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

size(block::Electron)
int64_t size() const;

dim

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

dim(block::Electron)
int64_tdim() const;

isreal

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

isreal(block::Electron; precision::Real=1e-12)
int64_t isreal(double precision = 1e-12) const;

Usage Example

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
p1 = Permutation([1, 2, 3, 4])
p2 = Permutation([2, 3, 4, 1])
p3 = Permutation([3, 4, 1, 2])
p4 = Permutation([4, 1, 2, 3])
group = PermutationGroup([p1, p2, p3, p4])
rep = Representation([1, -1, 1, -1])
block_sym = Electron(N, group, rep)
@show block_sym

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

@show n_sites(block_sym_np)
@show size(block_sym_np)

# Iteration
for pstate in block_sym_np
    @show pstate, index(block_sym_np, pstate)
end
@show permutation_group(block_sym_np)
@show irrep(block_sym_np)
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({1, -1, 1, -1});
auto block_sym = Electron(N, group, irrep);
XDIAG_SHOW(block_sym);

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

XDIAG_SHOW(block_sym_np.n_sites());
XDIAG_SHOW(block_sym_np.size());

// Iteration
for (auto pstate : block_sym_np) {
  Log("{} {}", to_string(pstate), block_sym_np.index(pstate));
}
XDIAG_SHOW(block_sym_np.permutation_group());
XDIAG_SHOW(block_sym_np.irrep());