TensorKit.jl aims to be a generic package for working with tensors as they appear throughout the physical sciences. TensorKit implements a parametric type Tensor (which is actually a specific case of the type TensorMap) and defines for these types a number of vector space operations (scalar multiplication, addition, norms and inner products), index operations (permutations) and linear algebra operations (multiplication, factorizations). Finally, tensor contractions can be performed using the @tensor macro from TensorOperations.jl.
Currently, most effort is oriented towards tensors as they appear in the context of quantum many body physics and in particular the field of tensor networks. Such tensors often have large dimensions and take on a specific structure when symmetries are present. To deal with generic symmetries, we employ notations and concepts from category theory all the way down to the definition of a tensor.
At the same time, TensorKit.jl focusses on computational efficiency and performance. The underlying storage of a tensor's data can be any DenseArray. Currently, certain operations are already multithreaded, either by distributing the different blocks in case of a structured tensor (i.e. with symmetries) or by using multithreading provided by the package Strided.jl. In the future, we also plan to investigate using CuArrays as underlying storage for the tensors data, so as to leverage GPUs for the different operations defined on tensors.
TensorKit.jl aims to be a generic package for working with tensors as they appear throughout the physical sciences. TensorKit implements a parametric type Tensor (which is actually a specific case of the type TensorMap) and defines for these types a number of vector space operations (scalar multiplication, addition, norms and inner products), index operations (permutations) and linear algebra operations (multiplication, factorizations). Finally, tensor contractions can be performed using the @tensor macro from TensorOperations.jl.
Currently, most effort is oriented towards tensors as they appear in the context of quantum many body physics and in particular the field of tensor networks. Such tensors often have large dimensions and take on a specific structure when symmetries are present. To deal with generic symmetries, we employ notations and concepts from category theory all the way down to the definition of a tensor.
At the same time, TensorKit.jl focusses on computational efficiency and performance. The underlying storage of a tensor's data can be any DenseArray. Currently, certain operations are already multithreaded, either by distributing the different blocks in case of a structured tensor (i.e. with symmetries) or by using multithreading provided by the package Strided.jl. In the future, we also plan to investigate using CuArrays as underlying storage for the tensors data, so as to leverage GPUs for the different operations defined on tensors.
Abstract type for representing the (isomorphism classes of) simple objects in (unitary and pivotal) (pre-)fusion categories, e.g. the irreducible representations of a finite or compact group. Subtypes I<:Sector as the set of labels of a GradedSpace.
Every new I<:Sector should implement the following methods:
one(::Type{I}): unit element of I
conj(a::I): $a̅$, conjugate or dual label of $a$
⊗(a::I, b::I): iterable with unique fusion outputs of $a ⊗ b$ (i.e. don't repeat in case of multiplicities)
Nsymbol(a::I, b::I, c::I): number of times c appears in a ⊗ b, i.e. the multiplicity
FusionStyle(::Type{I}): UniqueFusion(), SimpleFusion() or GenericFusion()
Singleton type to represent an iterator over the possible values of type I, whose instance is obtained as values(I). For a new I::Sector, the following should be defined
Base.iterate(::SectorValues{I}[, state]): iterate over the values
Base.IteratorSize(::Type{SectorValues{I}}): HasLenght(), SizeUnkown() or IsInfinite() depending on whether the number of values of type I is finite (and sufficiently small) or infinite; for a large number of values, SizeUnknown() is recommend because this will trigger the use of GenericGradedSpace.
If IteratorSize(I) == HasLength(), also the following must be implemented:
Base.length(::SectorValues{I}): the number of different values
Base.getindex(::SectorValues{I}, i::Int): a mapping between an index i and an instance of I
findindex(::SectorValues{I}, c::I): reverse mapping between a value c::I and an index i::Integer ∈ 1:length(values(I))
Return the type of fusion behavior of sectors of type I, which can be either
UniqueFusion(): single fusion output when fusing two sectors;
SimpleFusion(): multiple outputs, but every output occurs at most one, also known as multiplicity free (e.g. irreps of $SU(2)$);
GenericFusion(): multiple outputs that can occur more than once (e.g. irreps of $SU(3)$).
There is an abstract supertype MultipleFusion of which both SimpleFusion and GenericFusion are subtypes. Furthermore, there is a type alias MultiplicityFreeFusion for those fusion types which do not require muliplicity labels, i.e. MultiplicityFreeFusion = Union{UniqueFusion,SimpleFusion}.
Return the type of braiding and twist behavior of sectors of type I, which can be either
Bosonic(): symmetric braiding with trivial twist (i.e. identity)
Fermionic(): symmetric braiding with non-trivial twist (squares to identity)
Anyonic(): general $R_(a,b)^c$ phase or matrix (depending on SimpleFusion or GenericFusion fusion) and arbitrary twists
Note that Bosonic and Fermionic are subtypes of SymmetricBraiding, which means that braids are in fact equivalent to crossings (i.e. braiding twice is an identity: isone(Rsymbol(b,a,c)*Rsymbol(a,b,c)) == true) and permutations are uniquely defined.
abstract type AbstractIrrep{G<:Group} <: Sector end
Abstract supertype for sectors which corresponds to irreps (irreducible representations) of a group G. As we assume unitary representations, these would be finite groups or compact Lie groups. Note that this could also include projective rather than linear representations.
Actual concrete implementations of those irreps can be obtained as Irrep[G], or via their actual name, which generically takes the form (asciiG)Irrep, i.e. the ASCII spelling of the group name followed by Irrep.
All irreps have BraidingStyle equal to Bosonic() and thus trivial twists.
Represents irreps of the group $ℤ_N$ for some value of N<64. (We need 2*(N-1) <= 127 in order for a ⊗ b to work correctly.) For N equals 2, 3 or 4, ℤ{N} can be replaced by ℤ₂, ℤ₃, ℤ₄, whereas Parity is a synonym for Irrep{ℤ₂}. An arbitrary Integern can be provided to the constructor, but only the value mod(n, N) is relevant.
Represents irreps of the group $U₁$. The irrep is labelled by a charge, which should be an integer for a linear representation. However, it is often useful to allow half integers to represent irreps ofU₁subgroups ofSU₂, such as the Sz of spin-1/2 system. Hence, the charge is stored as aHalfIntfrom the package HalfIntegers.jl, but can be entered as arbitraryReal`. The sequence of the charges is: 0, 1/2, -1/2, 1, -1, ...
Represents irreps of the group $SU₂$. The irrep is labelled by a half integer j which can be entered as an abitrary Real, but is stored as a HalfInt from the HalfIntegers.jl package.
CU1Irrep(j, s = ifelse(j>zero(j), 2, 0))
-Irrep[CU₁](j, s = ifelse(j>zero(j), 2, 0))
Represents irreps of the group $U₁ ⋊ C$ ($U₁$ and charge conjugation or reflection), which is also known as just O₂. The irrep is labelled by a positive half integer j (the $U₁$ charge) and an integer s indicating the behaviour under charge conjugation. They take values:
if j == 0, s = 0 (trivial charge conjugation) or s = 1 (non-trivial charge conjugation)
Represents the anyons (isomorphism classes of simple objects) of the Fibonacci fusion category. It can take two values, corresponding to the trivial sector FibonacciAnyon(:I) and the non-trivial sector FibonacciAnyon(:τ) with fusion rules $τ ⊗ τ = 1 ⊕ τ$.
Represents a fusion tree of sectors of type I<:Sector, fusing (or splitting) N uncoupled sectors to a coupled sector. (It actually represents a splitting tree, but fusion tree is a more common term). The isdual field indicates whether an isomorphism is present (if the corresponding value is true) or not. The field uncoupled contains the sectors coming out of the splitting trees, before the possible 𝑍 isomorphism. This fusion tree has M=max(0, N-2) inner lines. Furthermore, for FusionStyle(I) isa GenericFusion, the L=max(0, N-1) corresponding vertices carry a label of type T. If FusionStyle(I) isa MultiplicityFreeFusion,T = Nothing`.
A constant of a singleton type used as Irrep[G] with G<:Group a type of group, to construct or obtain a concrete subtype of AbstractIrrep{G} that implements the data structure used to represent irreducible representations of the group G.
Nsymbol(a::I, b::I, c::I) where {I<:Sector} -> Integer
Return an Integer representing the number of times c appears in the fusion product a ⊗ b. Could be a Bool if FusionStyle(I) == UniqueFusion() or SimpleFusion().
Fsymbol(a::I, b::I, c::I, d::I, e::I, f::I) where {I<:Sector}
Return the F-symbol $F^{abc}_d$ that associates the two different fusion orders of sectors a, b and c into an ouput sector d, using either an intermediate sector $a ⊗ b → e$ or $b ⊗ c → f$:
a-<-μ-<-e-<-ν-<-d a-<-λ-<-d
+Symmetry sectors an fusion trees · TensorKit.jl
Abstract type for representing the (isomorphism classes of) simple objects in (unitary and pivotal) (pre-)fusion categories, e.g. the irreducible representations of a finite or compact group. Subtypes I<:Sector as the set of labels of a GradedSpace.
Every new I<:Sector should implement the following methods:
one(::Type{I}): unit element of I
conj(a::I): $a̅$, conjugate or dual label of $a$
⊗(a::I, b::I): iterable with unique fusion outputs of $a ⊗ b$ (i.e. don't repeat in case of multiplicities)
Nsymbol(a::I, b::I, c::I): number of times c appears in a ⊗ b, i.e. the multiplicity
FusionStyle(::Type{I}): UniqueFusion(), SimpleFusion() or GenericFusion()
Singleton type to represent an iterator over the possible values of type I, whose instance is obtained as values(I). For a new I::Sector, the following should be defined
Base.iterate(::SectorValues{I}[, state]): iterate over the values
Base.IteratorSize(::Type{SectorValues{I}}): HasLenght(), SizeUnkown() or IsInfinite() depending on whether the number of values of type I is finite (and sufficiently small) or infinite; for a large number of values, SizeUnknown() is recommend because this will trigger the use of GenericGradedSpace.
If IteratorSize(I) == HasLength(), also the following must be implemented:
Base.length(::SectorValues{I}): the number of different values
Base.getindex(::SectorValues{I}, i::Int): a mapping between an index i and an instance of I
findindex(::SectorValues{I}, c::I): reverse mapping between a value c::I and an index i::Integer ∈ 1:length(values(I))
Return the type of fusion behavior of sectors of type I, which can be either
UniqueFusion(): single fusion output when fusing two sectors;
SimpleFusion(): multiple outputs, but every output occurs at most one, also known as multiplicity free (e.g. irreps of $SU(2)$);
GenericFusion(): multiple outputs that can occur more than once (e.g. irreps of $SU(3)$).
There is an abstract supertype MultipleFusion of which both SimpleFusion and GenericFusion are subtypes. Furthermore, there is a type alias MultiplicityFreeFusion for those fusion types which do not require muliplicity labels, i.e. MultiplicityFreeFusion = Union{UniqueFusion,SimpleFusion}.
Return the type of braiding and twist behavior of sectors of type I, which can be either
Bosonic(): symmetric braiding with trivial twist (i.e. identity)
Fermionic(): symmetric braiding with non-trivial twist (squares to identity)
Anyonic(): general $R_(a,b)^c$ phase or matrix (depending on SimpleFusion or GenericFusion fusion) and arbitrary twists
Note that Bosonic and Fermionic are subtypes of SymmetricBraiding, which means that braids are in fact equivalent to crossings (i.e. braiding twice is an identity: isone(Rsymbol(b,a,c)*Rsymbol(a,b,c)) == true) and permutations are uniquely defined.
abstract type AbstractIrrep{G<:Group} <: Sector end
Abstract supertype for sectors which corresponds to irreps (irreducible representations) of a group G. As we assume unitary representations, these would be finite groups or compact Lie groups. Note that this could also include projective rather than linear representations.
Actual concrete implementations of those irreps can be obtained as Irrep[G], or via their actual name, which generically takes the form (asciiG)Irrep, i.e. the ASCII spelling of the group name followed by Irrep.
All irreps have BraidingStyle equal to Bosonic() and thus trivial twists.
Represents irreps of the group $ℤ_N$ for some value of N<64. (We need 2*(N-1) <= 127 in order for a ⊗ b to work correctly.) For N equals 2, 3 or 4, ℤ{N} can be replaced by ℤ₂, ℤ₃, ℤ₄, whereas Parity is a synonym for Irrep{ℤ₂}. An arbitrary Integern can be provided to the constructor, but only the value mod(n, N) is relevant.
Represents irreps of the group $U₁$. The irrep is labelled by a charge, which should be an integer for a linear representation. However, it is often useful to allow half integers to represent irreps ofU₁subgroups ofSU₂, such as the Sz of spin-1/2 system. Hence, the charge is stored as aHalfIntfrom the package HalfIntegers.jl, but can be entered as arbitraryReal`. The sequence of the charges is: 0, 1/2, -1/2, 1, -1, ...
Represents irreps of the group $SU₂$. The irrep is labelled by a half integer j which can be entered as an abitrary Real, but is stored as a HalfInt from the HalfIntegers.jl package.
CU1Irrep(j, s = ifelse(j>zero(j), 2, 0))
+Irrep[CU₁](j, s = ifelse(j>zero(j), 2, 0))
Represents irreps of the group $U₁ ⋊ C$ ($U₁$ and charge conjugation or reflection), which is also known as just O₂. The irrep is labelled by a positive half integer j (the $U₁$ charge) and an integer s indicating the behaviour under charge conjugation. They take values:
if j == 0, s = 0 (trivial charge conjugation) or s = 1 (non-trivial charge conjugation)
Represents the anyons (isomorphism classes of simple objects) of the Fibonacci fusion category. It can take two values, corresponding to the trivial sector FibonacciAnyon(:I) and the non-trivial sector FibonacciAnyon(:τ) with fusion rules $τ ⊗ τ = 1 ⊕ τ$.
Represents a fusion tree of sectors of type I<:Sector, fusing (or splitting) N uncoupled sectors to a coupled sector. (It actually represents a splitting tree, but fusion tree is a more common term). The isdual field indicates whether an isomorphism is present (if the corresponding value is true) or not. The field uncoupled contains the sectors coming out of the splitting trees, before the possible 𝑍 isomorphism. This fusion tree has M=max(0, N-2) inner lines. Furthermore, for FusionStyle(I) isa GenericFusion, the L=max(0, N-1) corresponding vertices carry a label of type T. If FusionStyle(I) isa MultiplicityFreeFusion,T = Nothing`.
A constant of a singleton type used as Irrep[G] with G<:Group a type of group, to construct or obtain a concrete subtype of AbstractIrrep{G} that implements the data structure used to represent irreducible representations of the group G.
Nsymbol(a::I, b::I, c::I) where {I<:Sector} -> Integer
Return an Integer representing the number of times c appears in the fusion product a ⊗ b. Could be a Bool if FusionStyle(I) == UniqueFusion() or SimpleFusion().
Fsymbol(a::I, b::I, c::I, d::I, e::I, f::I) where {I<:Sector}
Return the F-symbol $F^{abc}_d$ that associates the two different fusion orders of sectors a, b and c into an ouput sector d, using either an intermediate sector $a ⊗ b → e$ or $b ⊗ c → f$:
a-<-μ-<-e-<-ν-<-d a-<-λ-<-d
∨ ∨ -> Fsymbol(a,b,c,d,e,f)[μ,ν,κ,λ] ∨
b c f
v
b-<-κ
∨
- c
If FusionStyle(I) is UniqueFusion or SimpleFusion, the F-symbol is a number. Otherwise it is a rank 4 array of size (Nsymbol(a, b, e), Nsymbol(e, c, d), Nsymbol(b, c, f), Nsymbol(a, f, d)).
Returns the R-symbol $R^{ab}_c$ that maps between $c → a ⊗ b$ and $c → b ⊗ a$ as in
a -<-μ-<- c b -<-ν-<- c
+ c
If FusionStyle(I) is UniqueFusion or SimpleFusion, the F-symbol is a number. Otherwise it is a rank 4 array of size (Nsymbol(a, b, e), Nsymbol(e, c, d), Nsymbol(b, c, f), Nsymbol(a, f, d)).
Returns the R-symbol $R^{ab}_c$ that maps between $c → a ⊗ b$ and $c → b ⊗ a$ as in
a -<-μ-<- c b -<-ν-<- c
∨ -> Rsymbol(a,b,c)[μ,ν] v
- b a
If FusionStyle(I) is UniqueFusion() or SimpleFusion(), the R-symbol is a number. Otherwise it is a square matrix with row and column size Nsymbol(a,b,c) == Nsymbol(b,a,c).
Return the value of $B^{ab}_c$ which appears in transforming a splitting vertex into a fusion vertex using the transformation
a -<-μ-<- c a -<-ν-<- c
+ b a
If FusionStyle(I) is UniqueFusion() or SimpleFusion(), the R-symbol is a number. Otherwise it is a square matrix with row and column size Nsymbol(a,b,c) == Nsymbol(b,a,c).
Return the value of $B^{ab}_c$ which appears in transforming a splitting vertex into a fusion vertex using the transformation
a -<-μ-<- c a -<-ν-<- c
∨ -> √(dim(c)/dim(a)) * Bsymbol(a,b,c)[μ,ν] ∧
- b dual(b)
If FusionStyle(I) is UniqueFusion() or SimpleFusion(), the B-symbol is a number. Otherwise it is a square matrix with row and column size Nsymbol(a, b, c) == Nsymbol(c, dual(b), a).
vertex_ind2label(k::Int, a::I, b::I, c::I) where {I<:Sector}
Convert the index k of the fusion vertex (a,b)->c into a label. For FusionStyle(I) == UniqueFusion() or FusionStyle(I) == MultipleFusion(), where every fusion output occurs only once and k == 1, the default is to suppress vertex labels by setting them equal to nothing. For FusionStyle(I) == GenericFusion(), the default is to just use k, unless a specialized method is provided.
Given two sectors s₁ and s₂, which label an isomorphism class of simple objects in a fusion category $C₁$ and $C₂$, s1 ⊠ s2 (obtained as oxtimes+TAB) labels the isomorphism class of simple objects in the Deligne tensor product category $C₁ ⊠ C₂$.
The Deligne tensor product also works in the type domain and for spaces and tensors. For group representations, we have Irrep[G₁] ⊠ Irrep[G₂] == Irrep[G₁ × G₂].
braid(f₁::FusionTree{I}, f₂::FusionTree{I},
+ b dual(b)
If FusionStyle(I) is UniqueFusion() or SimpleFusion(), the B-symbol is a number. Otherwise it is a square matrix with row and column size Nsymbol(a, b, c) == Nsymbol(c, dual(b), a).
vertex_ind2label(k::Int, a::I, b::I, c::I) where {I<:Sector}
Convert the index k of the fusion vertex (a,b)->c into a label. For FusionStyle(I) == UniqueFusion() or FusionStyle(I) == MultipleFusion(), where every fusion output occurs only once and k == 1, the default is to suppress vertex labels by setting them equal to nothing. For FusionStyle(I) == GenericFusion(), the default is to just use k, unless a specialized method is provided.
Given two sectors s₁ and s₂, which label an isomorphism class of simple objects in a fusion category $C₁$ and $C₂$, s1 ⊠ s2 (obtained as oxtimes+TAB) labels the isomorphism class of simple objects in the Deligne tensor product category $C₁ ⊠ C₂$.
The Deligne tensor product also works in the type domain and for spaces and tensors. For group representations, we have Irrep[G₁] ⊠ Irrep[G₂] == Irrep[G₁ × G₂].
Input is a fusion-splitting tree pair that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the splitting tree f₁ and fusion tree f₂, such that the incoming sectors f₂.uncoupled are fused to f₁.coupled == f₂.coupled and then to the outgoing sectors f₁.uncoupled. Compute new trees and corresponding coefficients obtained from repartitioning and braiding the tree such that sectors p1 become outgoing and sectors p2 become incoming. The uncoupled indices in splitting tree f₁ and fusion tree f₂ have levels (or depths) levels1 and levels2 respectively, which determines how indices braid. In particular, if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Input is a fusion-splitting tree pair that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the splitting tree f₁ and fusion tree f₂, such that the incoming sectors f₂.uncoupled are fused to f₁.coupled == f₂.coupled and then to the outgoing sectors f₁.uncoupled. Compute new trees and corresponding coefficients obtained from repartitioning and braiding the tree such that sectors p1 become outgoing and sectors p2 become incoming. The uncoupled indices in splitting tree f₁ and fusion tree f₂ have levels (or depths) levels1 and levels2 respectively, which determines how indices braid. In particular, if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Input is a double fusion tree that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the individual trees of outgoing (t1) and incoming sectors (t2) respectively (with identical coupled sector t1.coupled == t2.coupled). Computes new trees and corresponding coefficients obtained from repartitioning and permuting the tree such that sectors p1 become outgoing and sectors p2 become incoming.
Perform a braiding of the uncoupled indices of the fusion tree f and return the result as a <:AbstractDict of output trees and corresponding coefficients. The braiding is determined by specifying that the new sector at position k corresponds to the sector that was originally at the position i = p[k], and assigning to every index i of the original fusion tree a distinct level or depth levels[i]. This permutation is then decomposed into elementary swaps between neighbouring indices, where the swaps are applied as braids such that if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Perform a permutation of the uncoupled indices of the fusion tree f and returns the result as a <:AbstractDict of output trees and corresponding coefficients; this requires that BraidingStyle(sectortype(f)) isa SymmetricBraiding.
Input is a double fusion tree that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the individual trees of outgoing (f₁) and incoming sectors (f₂) respectively (with identical coupled sector f₁.coupled == f₂.coupled). Computes new trees and corresponding coefficients obtained from repartitioning the tree by bending incoming to outgoing sectors (or vice versa) in order to have N outgoing sectors.
Perform an elementary braid (Artin generator) of neighbouring uncoupled indices i and i+1 on a fusion tree f, and returns the result as a dictionary of output trees and corresponding coefficients.
The keyword inv determines whether index i will braid above or below index i+1, i.e. applying artin_braid(f′, i; inv = true) to all the outputs f′ of artin_braid(f, i; inv = false) and collecting the results should yield a single fusion tree with non-zero coefficient, namely f with coefficient 1. This keyword has no effect if BraidingStyle(sectortype(f)) isa SymmetricBraiding.
Attach a fusion tree f₂ to the uncoupled leg i of the fusion tree f₁ and bring it into a linear combination of fusion trees in standard form. This requires that f₂.coupled == f₁.uncoupled[i] and f₁.isdual[i] == false.
Split a fusion tree into two. The first tree has as uncoupled sectors the first M uncoupled sectors of the input tree f, whereas its coupled sector corresponds to the internal sector between uncoupled sectors M and M+1 of the original tree f. The second tree has as first uncoupled sector that same internal sector of f, followed by remaining N-M uncoupled sectors of f. It couples to the same sector as f. This operation is the inverse of insertat in the sense that if f₁, f₂ = split(t, M) ⇒ f == insertat(f₂, 1, f₁).
Merge two fusion trees together to a linear combination of fusion trees whose uncoupled sectors are those of f₁ followed by those of f₂, and where the two coupled sectors of f₁ and f₂ are further fused to c. In case of FusionStyle(I) == GenericFusion(), also a degeneracy label μ for the fusion of the coupled sectors of f₁ and f₂ to c needs to be specified.
Input is a double fusion tree that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the individual trees of outgoing (t1) and incoming sectors (t2) respectively (with identical coupled sector t1.coupled == t2.coupled). Computes new trees and corresponding coefficients obtained from repartitioning and permuting the tree such that sectors p1 become outgoing and sectors p2 become incoming.
Perform a braiding of the uncoupled indices of the fusion tree f and return the result as a <:AbstractDict of output trees and corresponding coefficients. The braiding is determined by specifying that the new sector at position k corresponds to the sector that was originally at the position i = p[k], and assigning to every index i of the original fusion tree a distinct level or depth levels[i]. This permutation is then decomposed into elementary swaps between neighbouring indices, where the swaps are applied as braids such that if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Perform a permutation of the uncoupled indices of the fusion tree f and returns the result as a <:AbstractDict of output trees and corresponding coefficients; this requires that BraidingStyle(sectortype(f)) isa SymmetricBraiding.
Input is a double fusion tree that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the individual trees of outgoing (f₁) and incoming sectors (f₂) respectively (with identical coupled sector f₁.coupled == f₂.coupled). Computes new trees and corresponding coefficients obtained from repartitioning the tree by bending incoming to outgoing sectors (or vice versa) in order to have N outgoing sectors.
Perform an elementary braid (Artin generator) of neighbouring uncoupled indices i and i+1 on a fusion tree f, and returns the result as a dictionary of output trees and corresponding coefficients.
The keyword inv determines whether index i will braid above or below index i+1, i.e. applying artin_braid(f′, i; inv = true) to all the outputs f′ of artin_braid(f, i; inv = false) and collecting the results should yield a single fusion tree with non-zero coefficient, namely f with coefficient 1. This keyword has no effect if BraidingStyle(sectortype(f)) isa SymmetricBraiding.
Attach a fusion tree f₂ to the uncoupled leg i of the fusion tree f₁ and bring it into a linear combination of fusion trees in standard form. This requires that f₂.coupled == f₁.uncoupled[i] and f₁.isdual[i] == false.
Split a fusion tree into two. The first tree has as uncoupled sectors the first M uncoupled sectors of the input tree f, whereas its coupled sector corresponds to the internal sector between uncoupled sectors M and M+1 of the original tree f. The second tree has as first uncoupled sector that same internal sector of f, followed by remaining N-M uncoupled sectors of f. It couples to the same sector as f. This operation is the inverse of insertat in the sense that if f₁, f₂ = split(t, M) ⇒ f == insertat(f₂, 1, f₁).
Merge two fusion trees together to a linear combination of fusion trees whose uncoupled sectors are those of f₁ followed by those of f₂, and where the two coupled sectors of f₁ and f₂ are further fused to c. In case of FusionStyle(I) == GenericFusion(), also a degeneracy label μ for the fusion of the coupled sectors of f₁ and f₂ to c needs to be specified.
Abstract type at the top of the type hierarchy for denoting fields over which vector spaces (or more generally, linear categories) can be defined. Two common fields are ℝ and ℂ, representing the field of real or complex numbers respectively.
Abstract type at the top of the type hierarchy for denoting vector spaces, or, more accurately, 𝕜-linear categories. All instances of subtypes of VectorSpace will represent objects in 𝕜-linear monoidal categories.
Elementary finite-dimensional vector space over a field that can be used as the index space corresponding to the indices of a tensor. ElementarySpace is a supertype for all vector spaces (objects) that can be associated with the individual indices of a tensor, as hinted to by its alias IndexSpace.
Every elementary vector space should respond to the methods conj and dual, returning the complex conjugate space and the dual space respectively. The complex conjugate of the dual space is obtained as dual(conj(V)) === conj(dual(V)). These different spaces should be of the same type, so that a tensor can be defined as an element of a homogeneous tensor product of these spaces.
A finite-dimensional space over an arbitrary field 𝕜 without additional structure. It is thus characterized by its dimension, and whether or not it is the dual and/or conjugate space. For a real field 𝕜, the space and its conjugate are the same.
A real Euclidean space ℝ^d, which is therefore self-dual. CartesianSpace has no additonal structure and is completely characterised by its dimension d. This is the vector space that is implicitly assumed in most of matrix algebra.
A standard complex vector space ℂ^d with Euclidean inner product and no additional structure. It is completely characterised by its dimension and whether its the normal space or its dual (which is canonically isomorphic to the conjugate space).
Abstract type at the top of the type hierarchy for denoting fields over which vector spaces (or more generally, linear categories) can be defined. Two common fields are ℝ and ℂ, representing the field of real or complex numbers respectively.
Abstract type at the top of the type hierarchy for denoting vector spaces, or, more accurately, 𝕜-linear categories. All instances of subtypes of VectorSpace will represent objects in 𝕜-linear monoidal categories.
Elementary finite-dimensional vector space over a field that can be used as the index space corresponding to the indices of a tensor. ElementarySpace is a supertype for all vector spaces (objects) that can be associated with the individual indices of a tensor, as hinted to by its alias IndexSpace.
Every elementary vector space should respond to the methods conj and dual, returning the complex conjugate space and the dual space respectively. The complex conjugate of the dual space is obtained as dual(conj(V)) === conj(dual(V)). These different spaces should be of the same type, so that a tensor can be defined as an element of a homogeneous tensor product of these spaces.
A finite-dimensional space over an arbitrary field 𝕜 without additional structure. It is thus characterized by its dimension, and whether or not it is the dual and/or conjugate space. For a real field 𝕜, the space and its conjugate are the same.
A real Euclidean space ℝ^d, which is therefore self-dual. CartesianSpace has no additonal structure and is completely characterised by its dimension d. This is the vector space that is implicitly assumed in most of matrix algebra.
A standard complex vector space ℂ^d with Euclidean inner product and no additional structure. It is completely characterised by its dimension and whether its the normal space or its dual (which is canonically isomorphic to the conjugate space).
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A ProductSpace is a tensor product space of N vector spaces of type S<:ElementarySpace. Only tensor products between ElementarySpace objects of the same type are allowed.
A constant of a singleton type used as Vect[I] with I<:Sector a type of sector, to construct or obtain the concrete type GradedSpace{I,D} instances without having to specify D.
A constant of a singleton type used as Rep[G] with G<:Group a type of group, to construct or obtain the concrete type GradedSpace{Irrep[G],D} instances without having to specify D. Note that Rep[G] == Vect[Irrep[G]].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A ProductSpace is a tensor product space of N vector spaces of type S<:ElementarySpace. Only tensor products between ElementarySpace objects of the same type are allowed.
A constant of a singleton type used as Vect[I] with I<:Sector a type of sector, to construct or obtain the concrete type GradedSpace{I,D} instances without having to specify D.
A constant of a singleton type used as Rep[G] with G<:Group a type of group, to construct or obtain the concrete type GradedSpace{Irrep[G],D} instances without having to specify D. Note that Rep[G] == Vect[Irrep[G]].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
sectors(P::ProductSpace{S, N}) where {S<:ElementarySpace}
Return an iterator over all possible combinations of sectors (represented as an NTuple{N, sectortype(S)}) that can appear within the tensor product space P.
Return an iterator over the different unique coupled sector labels, i.e. the different fusion outputs that can be obtained by fusing the sectors present in the different spaces that make up the ProductSpace instance.
Return an iterator over the different unique coupled sector labels, i.e. the intersection of the different fusion outputs that can be obtained by fusing the sectors present in the domain, as well as from the codomain.
Return the total dimension of a coupled sector c in the product space, by summing over all dim(P, s) for all tuples of sectors s::NTuple{N, <:Sector} that can fuse to c, counted with the correct multiplicity (i.e. number of ways in which s can fuse to c).
Return a single vector space of type S that has the same value of isdual as dual(V), but yet is isomorphic to V rather than to dual(V). The spaces flip(V) and dual(V) only differ in the case of GradedSpace{I}.
⊕(V₁::S, V₂::S, V₃::S...) where {S<:ElementarySpace} -> S
Return the corresponding vector space of type S that represents the direct sum sum of the spaces V₁, V₂, ... Note that all the individual spaces should have the same value for isdual, as otherwise the direct sum is not defined.
Return the corresponding vector space of type S that represents the trivial one-dimensional space, i.e. the space that is isomorphic to the corresponding field. Note that this is different from one(V::S), which returns the empty product space ProductSpace{S,0}(()).
Return the supremum of a number of elementary spaces, i.e. an instance V::ElementarySpace such that V ≿ V₁, V ≿ V₂, ... and no other W ≺ V has this property. This requires that all arguments have the same value of isdual( ), and also the return value V will have the same value.
Return the infimum of a number of elementary spaces, i.e. an instance V::ElementarySpace such that V ≾ V₁, V ≾ V₂, ... and no other W ≻ V has this property. This requires that all arguments have the same value of isdual( ), and also the return value V will have the same value.
⊗(V₁::S, V₂::S, V₃::S...) where {S<:ElementarySpace} -> S
Create a ProductSpace{S}(V₁, V₂, V₃...) representing the tensor product of several elementary vector spaces. For convience, Julia's regular multiplication operator * applied to vector spaces has the same effect.
The tensor product structure is preserved, see fuse for returning a single elementary space of type S that is isomorphic to this tensor product.
Given two sectors s₁ and s₂, which label an isomorphism class of simple objects in a fusion category $C₁$ and $C₂$, s1 ⊠ s2 (obtained as oxtimes+TAB) labels the isomorphism class of simple objects in the Deligne tensor product category $C₁ ⊠ C₂$.
The Deligne tensor product also works in the type domain and for spaces and tensors. For group representations, we have Irrep[G₁] ⊠ Irrep[G₂] == Irrep[G₁ × G₂].
Given two vector spaces V₁ and V₂ (ElementarySpace or ProductSpace), or thus, objects of corresponding fusion categories $C₁$ and $C₂$, $V₁ ⊠ V₂$ constructs the Deligne tensor product, an object in $C₁ ⊠ C₂$ which is the natural tensor product of those categories. In particular, the corresponding type of sectors (simple objects) is given by sectortype(V₁ ⊠ V₂) == sectortype(V₁) ⊠ sectortype(V₂) and can be thought of as a tuple of the individual sectors.
The Deligne tensor product also works in the type domain and for sectors and tensors. For group representations, we have Rep[G₁] ⊠ Rep[G₂] == Rep[G₁ × G₂], i.e. these are the natural representation spaces of the direct product of two groups.
For P::ProductSpace{S,N}, this adds an extra tensor product factor at position 1 <= i <= N+1 (last position by default) which is just the S-equivalent of the underlying field of scalars, i.e. oneunit(S). With the keyword arguments, one can choose to insert the conjugated or dual space instead, which are all isomorphic to the field of scalars.
This document was generated with Documenter.jl version 0.27.25 on Wednesday 1 November 2023. Using Julia version 1.9.3.
+end
A complex Euclidean space with a direct sum structure corresponding to labels in a set I, the objects of which have the structure of a monoid with respect to a monoidal product ⊗. In practice, we restrict the label set to be a set of superselection sectors of type I<:Sector, e.g. the set of distinct irreps of a finite or compact group, or the isomorphism classes of simple objects of a unitary and pivotal (pre-)fusion category.
Here dims represents the degeneracy or multiplicity of every sector.
The data structure D of dims will depend on the result Base.IteratorElsize(values(I)); if the result is of type HasLength or HasShape, dims will be stored in a NTuple{N,Int} with N = length(values(I)). This requires that a sector s::I can be transformed into an index via s == getindex(values(I), i) and i == findindex(values(I), s). If Base.IteratorElsize(values(I)) results IsInfinite() or SizeUnknown(), a SectorDict{I,Int} is used to store the non-zero degeneracy dimensions with the corresponding sector as key. The parameter D is hidden from the user and should typically be of no concern.
The concrete type GradedSpace{I,D} with correct D can be obtained as Vect[I], or if I == Irrep[G] for some G<:Group, as Rep[G].
sectors(P::ProductSpace{S, N}) where {S<:ElementarySpace}
Return an iterator over all possible combinations of sectors (represented as an NTuple{N, sectortype(S)}) that can appear within the tensor product space P.
Return an iterator over the different unique coupled sector labels, i.e. the different fusion outputs that can be obtained by fusing the sectors present in the different spaces that make up the ProductSpace instance.
Return an iterator over the different unique coupled sector labels, i.e. the intersection of the different fusion outputs that can be obtained by fusing the sectors present in the domain, as well as from the codomain.
Return the total dimension of a coupled sector c in the product space, by summing over all dim(P, s) for all tuples of sectors s::NTuple{N, <:Sector} that can fuse to c, counted with the correct multiplicity (i.e. number of ways in which s can fuse to c).
Return a single vector space of type S that has the same value of isdual as dual(V), but yet is isomorphic to V rather than to dual(V). The spaces flip(V) and dual(V) only differ in the case of GradedSpace{I}.
⊕(V₁::S, V₂::S, V₃::S...) where {S<:ElementarySpace} -> S
Return the corresponding vector space of type S that represents the direct sum sum of the spaces V₁, V₂, ... Note that all the individual spaces should have the same value for isdual, as otherwise the direct sum is not defined.
Return the corresponding vector space of type S that represents the trivial one-dimensional space, i.e. the space that is isomorphic to the corresponding field. Note that this is different from one(V::S), which returns the empty product space ProductSpace{S,0}(()).
Return the supremum of a number of elementary spaces, i.e. an instance V::ElementarySpace such that V ≿ V₁, V ≿ V₂, ... and no other W ≺ V has this property. This requires that all arguments have the same value of isdual( ), and also the return value V will have the same value.
Return the infimum of a number of elementary spaces, i.e. an instance V::ElementarySpace such that V ≾ V₁, V ≾ V₂, ... and no other W ≻ V has this property. This requires that all arguments have the same value of isdual( ), and also the return value V will have the same value.
⊗(V₁::S, V₂::S, V₃::S...) where {S<:ElementarySpace} -> S
Create a ProductSpace{S}(V₁, V₂, V₃...) representing the tensor product of several elementary vector spaces. For convience, Julia's regular multiplication operator * applied to vector spaces has the same effect.
The tensor product structure is preserved, see fuse for returning a single elementary space of type S that is isomorphic to this tensor product.
Given two sectors s₁ and s₂, which label an isomorphism class of simple objects in a fusion category $C₁$ and $C₂$, s1 ⊠ s2 (obtained as oxtimes+TAB) labels the isomorphism class of simple objects in the Deligne tensor product category $C₁ ⊠ C₂$.
The Deligne tensor product also works in the type domain and for spaces and tensors. For group representations, we have Irrep[G₁] ⊠ Irrep[G₂] == Irrep[G₁ × G₂].
Given two vector spaces V₁ and V₂ (ElementarySpace or ProductSpace), or thus, objects of corresponding fusion categories $C₁$ and $C₂$, $V₁ ⊠ V₂$ constructs the Deligne tensor product, an object in $C₁ ⊠ C₂$ which is the natural tensor product of those categories. In particular, the corresponding type of sectors (simple objects) is given by sectortype(V₁ ⊠ V₂) == sectortype(V₁) ⊠ sectortype(V₂) and can be thought of as a tuple of the individual sectors.
The Deligne tensor product also works in the type domain and for sectors and tensors. For group representations, we have Rep[G₁] ⊠ Rep[G₂] == Rep[G₁ × G₂], i.e. these are the natural representation spaces of the direct product of two groups.
For P::ProductSpace{S,N}, this adds an extra tensor product factor at position 1 <= i <= N+1 (last position by default) which is just the S-equivalent of the underlying field of scalars, i.e. oneunit(S). With the keyword arguments, one can choose to insert the conjugated or dual space instead, which are all isomorphic to the field of scalars.
abstract type AbstractTensorMap{S<:IndexSpace, N₁, N₂} end
Abstract supertype of all tensor maps, i.e. linear maps between tensor products of vector spaces of type S<:IndexSpace. An AbstractTensorMap maps from an input space of type ProductSpace{S, N₂} to an output space of type ProductSpace{S, N₁}.
AbstractTensor{S<:IndexSpace, N} = AbstractTensorMap{S, N, 0}
Abstract supertype of all tensors, i.e. elements in the tensor product space of type ProductSpace{S, N}, built from elementary spaces of type S<:IndexSpace.
An AbstractTensor{S, N} is actually a special case AbstractTensorMap{S, N, 0}, i.e. a tensor map with only a non-trivial output space.
Specific subtype of AbstractTensorMap for representing tensor maps (morphisms in a tensor category) whose data is stored in blocks of some subtype of DenseMatrix.
Construct the identity endomorphism on space space, i.e. return a t::TensorMap with domain(t) == codomain(t) == V, where storagetype(t) = A can be specified.
abstract type AbstractTensorMap{S<:IndexSpace, N₁, N₂} end
Abstract supertype of all tensor maps, i.e. linear maps between tensor products of vector spaces of type S<:IndexSpace. An AbstractTensorMap maps from an input space of type ProductSpace{S, N₂} to an output space of type ProductSpace{S, N₁}.
AbstractTensor{S<:IndexSpace, N} = AbstractTensorMap{S, N, 0}
Abstract supertype of all tensors, i.e. elements in the tensor product space of type ProductSpace{S, N}, built from elementary spaces of type S<:IndexSpace.
An AbstractTensor{S, N} is actually a special case AbstractTensorMap{S, N, 0}, i.e. a tensor map with only a non-trivial output space.
Specific subtype of AbstractTensorMap for representing tensor maps (morphisms in a tensor category) whose data is stored in blocks of some subtype of DenseMatrix.
Construct the identity endomorphism on space space, i.e. return a t::TensorMap with domain(t) == codomain(t) == V, where storagetype(t) = A can be specified.
Return a t::TensorMap that implements a specific isomorphism between the codomain cod and the domain dom, and for which storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that isomorphism(cod, dom) == inv(isomorphism(dom, cod)).
See also unitary when InnerProductStyle(cod) === EuclideanProduct().
Return a t::TensorMap that implements a specific unitary isomorphism between the codomain cod and the domain dom, for which spacetype(dom) (== spacetype(cod)) must have an inner product. Furthermore, storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that unitary(cod, dom) == inv(unitary(dom, cod)) = adjoint(unitary(dom, cod)).
Return a t::TensorMap that implements a specific isometry that embeds the domain dom into the codomain cod, and which requires that spacetype(dom) (== spacetype(cod)) has an Euclidean inner product. An isometry t is such that its adjoint t' is the left inverse of t, i.e. t'*t = id(dom), while t*t' is some idempotent endomorphism of cod, i.e. it squares to itself. When dom and cod do not allow for such an isometric inclusion, an error will be thrown.
Return a t::TensorMap that implements a specific isomorphism between the codomain cod and the domain dom, and for which storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that isomorphism(cod, dom) == inv(isomorphism(dom, cod)).
See also unitary when InnerProductStyle(cod) === EuclideanProduct().
Return a t::TensorMap that implements a specific unitary isomorphism between the codomain cod and the domain dom, for which spacetype(dom) (== spacetype(cod)) must have an inner product. Furthermore, storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that unitary(cod, dom) == inv(unitary(dom, cod)) = adjoint(unitary(dom, cod)).
Return a t::TensorMap that implements a specific isometry that embeds the domain dom into the codomain cod, and which requires that spacetype(dom) (== spacetype(cod)) has an Euclidean inner product. An isometry t is such that its adjoint t' is the left inverse of t, i.e. t'*t = id(dom), while t*t' is some idempotent endomorphism of cod, i.e. it squares to itself. When dom and cod do not allow for such an isometric inclusion, an error will be thrown.
permute!(tdst::AbstractTensorMap{S,N₁,N₂}, tsrc::AbstractTensorMap{S},
(p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}) where {S,N₁,N₂}
- -> tdst
Write into tdst the result of permuting the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively.
See permute for creating a new tensor and add_permute! for a more general version.
Perform a braiding of the uncoupled indices of the fusion tree f and return the result as a <:AbstractDict of output trees and corresponding coefficients. The braiding is determined by specifying that the new sector at position k corresponds to the sector that was originally at the position i = p[k], and assigning to every index i of the original fusion tree a distinct level or depth levels[i]. This permutation is then decomposed into elementary swaps between neighbouring indices, where the swaps are applied as braids such that if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Write into tdst the result of permuting the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively.
See permute for creating a new tensor and add_permute! for a more general version.
Perform a braiding of the uncoupled indices of the fusion tree f and return the result as a <:AbstractDict of output trees and corresponding coefficients. The braiding is determined by specifying that the new sector at position k corresponds to the sector that was originally at the position i = p[k], and assigning to every index i of the original fusion tree a distinct level or depth levels[i]. This permutation is then decomposed into elementary swaps between neighbouring indices, where the swaps are applied as braids such that if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Input is a fusion-splitting tree pair that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the splitting tree f₁ and fusion tree f₂, such that the incoming sectors f₂.uncoupled are fused to f₁.coupled == f₂.coupled and then to the outgoing sectors f₁.uncoupled. Compute new trees and corresponding coefficients obtained from repartitioning and braiding the tree such that sectors p1 become outgoing and sectors p2 become incoming. The uncoupled indices in splitting tree f₁ and fusion tree f₂ have levels (or depths) levels1 and levels2 respectively, which determines how indices braid. In particular, if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Input is a fusion-splitting tree pair that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the splitting tree f₁ and fusion tree f₂, such that the incoming sectors f₂.uncoupled are fused to f₁.coupled == f₂.coupled and then to the outgoing sectors f₁.uncoupled. Compute new trees and corresponding coefficients obtained from repartitioning and braiding the tree such that sectors p1 become outgoing and sectors p2 become incoming. The uncoupled indices in splitting tree f₁ and fusion tree f₂ have levels (or depths) levels1 and levels2 respectively, which determines how indices braid. In particular, if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.
Return tensor tdst obtained by braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.
If copy=false, tdst might share data with tsrc whenever possible. Otherwise, a copy is always made.
To braid into an existing destination, see braid! and add_braid!
Return tensor tdst obtained by braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.
If copy=false, tdst might share data with tsrc whenever possible. Otherwise, a copy is always made.
To braid into an existing destination, see braid! and add_braid!
braid!(tdst::AbstractTensorMap{S,N₁,N₂}, tsrc::AbstractTensorMap{S},
(p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}, levels::Tuple) where {S,N₁,N₂}
- -> tdst
Write into tdst the result of braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.
See braid for creating a new tensor and add_braid! for a more general version.
leftorth(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple;
- alg::OrthogonalFactorizationAlgorithm = QRpos()) -> Q, R
Create orthonormal basis Q for indices in leftind, and remainder R such that permute(t, (leftind, rightind)) = Q*R.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftorth!(t, alg = QRpos()).
Different algorithms are available, namely QR(), QRpos(), SVD() and Polar(). QR() and QRpos() use a standard QR decomposition, producing an upper triangular matrix R. Polar() produces a Hermitian and positive semidefinite R. QRpos() corrects the standard QR decomposition such that the diagonal elements of R are positive. Only QRpos() and Polar() are unique (no residual freedom) so that they always return the same result for the same input tensor t.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and leftorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
Create orthonormal basis Q for indices in rightind, and remainder L such that permute(t, (leftind, rightind)) = L*Q.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightorth!(t, alg = LQpos()).
Different algorithms are available, namely LQ(), LQpos(), RQ(), RQpos(), SVD() and Polar(). LQ() and LQpos() produce a lower triangular matrix L and are computed using a QR decomposition of the transpose. RQ() and RQpos() produce an upper triangular remainder L and only works if the total left dimension is smaller than or equal to the total right dimension. LQpos() and RQpos() add an additional correction such that the diagonal elements of L are positive. Polar() produces a Hermitian and positive semidefinite L. Only LQpos(), RQpos() and Polar() are unique (no residual freedom) so that they always return the same result for the same input tensor t.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and rightorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
leftnull(t::AbstractTensor, (leftind, rightind)::Index2Tuple;
- alg::OrthogonalFactorizationAlgorithm = QRpos()) -> N
Create orthonormal basis for the orthogonal complement of the support of the indices in leftind, such that N' * permute(t, (leftind, rightind)) = 0.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftnull!(t, alg = QRpos()).
Different algorithms are available, namely QR() (or equivalently, QRpos()), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and leftnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
Write into tdst the result of braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.
See braid for creating a new tensor and add_braid! for a more general version.
leftorth(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple;
+ alg::OrthogonalFactorizationAlgorithm = QRpos()) -> Q, R
Create orthonormal basis Q for indices in leftind, and remainder R such that permute(t, (leftind, rightind)) = Q*R.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftorth!(t, alg = QRpos()).
Different algorithms are available, namely QR(), QRpos(), SVD() and Polar(). QR() and QRpos() use a standard QR decomposition, producing an upper triangular matrix R. Polar() produces a Hermitian and positive semidefinite R. QRpos() corrects the standard QR decomposition such that the diagonal elements of R are positive. Only QRpos() and Polar() are unique (no residual freedom) so that they always return the same result for the same input tensor t.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and leftorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
Create orthonormal basis Q for indices in rightind, and remainder L such that permute(t, (leftind, rightind)) = L*Q.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightorth!(t, alg = LQpos()).
Different algorithms are available, namely LQ(), LQpos(), RQ(), RQpos(), SVD() and Polar(). LQ() and LQpos() produce a lower triangular matrix L and are computed using a QR decomposition of the transpose. RQ() and RQpos() produce an upper triangular remainder L and only works if the total left dimension is smaller than or equal to the total right dimension. LQpos() and RQpos() add an additional correction such that the diagonal elements of L are positive. Polar() produces a Hermitian and positive semidefinite L. Only LQpos(), RQpos() and Polar() are unique (no residual freedom) so that they always return the same result for the same input tensor t.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and rightorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
leftnull(t::AbstractTensor, (leftind, rightind)::Index2Tuple;
+ alg::OrthogonalFactorizationAlgorithm = QRpos()) -> N
Create orthonormal basis for the orthogonal complement of the support of the indices in leftind, such that N' * permute(t, (leftind, rightind)) = 0.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftnull!(t, alg = QRpos()).
Different algorithms are available, namely QR() (or equivalently, QRpos()), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and leftnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
Create orthonormal basis for the orthogonal complement of the support of the indices in rightind, such that permute(t, (leftind, rightind))*N' = 0.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightnull!(t, alg = LQpos()).
Different algorithms are available, namely LQ() (or equivalently, LQpos), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and rightnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
tsvd(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple;
+ rtol::Real = eps(real(float(one(scalartype(t)))))*iszero(atol)) -> N
Create orthonormal basis for the orthogonal complement of the support of the indices in rightind, such that permute(t, (leftind, rightind))*N' = 0.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightnull!(t, alg = LQpos()).
Different algorithms are available, namely LQ() (or equivalently, LQpos), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and rightnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
Compute the (possibly truncated) singular value decomposition such that norm(permute(t, (leftind, rightind)) - U * S * V) ≈ ϵ, where ϵ thus represents the truncation error.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using tsvd!(t, trunc = notrunc(), p = 2).
A truncation parameter trunc can be specified for the new internal dimension, in which case a truncated singular value decomposition will be computed. Choices are:
notrunc(): no truncation (default);
truncerr(η::Real): truncates such that the p-norm of the truncated singular values is smaller than η times the p-norm of all singular values;
truncdim(χ::Int): truncates such that the equivalent total dimension of the internal vector space is no larger than χ;
truncspace(V): truncates such that the dimension of the internal vector space is smaller than that of V in any sector.
truncbelow(χ::Real): truncates such that every singular value is larger then χ ;
The method tsvd also returns the truncation error ϵ, computed as the p norm of the singular values that were truncated.
THe keyword alg can be equal to SVD() or SDD(), corresponding to the underlying LAPACK algorithm that computes the decomposition (_gesvd or _gesdd).
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and tsvd(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
eigh(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eigh assumes that the linear map is hermitian and D and V tensors with the same scalartype as t. See eig and eigen for non-hermitian tensors. Hermiticity requires that the tensor acts on inner product spaces, and the current implementation requires InnerProductStyle(t) === EuclideanProduct().
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigh!(t). Note that the permuted tensor on which eigh! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
eig(t::AbstractTensor, (leftind, rightind)::Index2Tuple; kwargs...) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eig assumes that the linear map is not hermitian and returns type stable complex valued D and V tensors for both real and complex valued t. See eigh for hermitian linear maps
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eig!(t). Note that the permuted tensor on which eig! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
permute(t, (leftind, rightind)) * V = V * D
Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.
eigen(t::AbstractTensor, (leftind, rightind)::Index2Tuple; kwargs...) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigen!(t). Note that the permuted tensor on which eigen! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
permute(t, (leftind, rightind)) * V = V * D
Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.
This document was generated with Documenter.jl version 0.27.25 on Wednesday 1 November 2023. Using Julia version 1.9.3.
+ -> U, S, V, ϵ
Compute the (possibly truncated) singular value decomposition such that norm(permute(t, (leftind, rightind)) - U * S * V) ≈ ϵ, where ϵ thus represents the truncation error.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using tsvd!(t, trunc = notrunc(), p = 2).
A truncation parameter trunc can be specified for the new internal dimension, in which case a truncated singular value decomposition will be computed. Choices are:
notrunc(): no truncation (default);
truncerr(η::Real): truncates such that the p-norm of the truncated singular values is smaller than η times the p-norm of all singular values;
truncdim(χ::Int): truncates such that the equivalent total dimension of the internal vector space is no larger than χ;
truncspace(V): truncates such that the dimension of the internal vector space is smaller than that of V in any sector.
truncbelow(χ::Real): truncates such that every singular value is larger then χ ;
The method tsvd also returns the truncation error ϵ, computed as the p norm of the singular values that were truncated.
THe keyword alg can be equal to SVD() or SDD(), corresponding to the underlying LAPACK algorithm that computes the decomposition (_gesvd or _gesdd).
Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and tsvd(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().
eigh(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eigh assumes that the linear map is hermitian and D and V tensors with the same scalartype as t. See eig and eigen for non-hermitian tensors. Hermiticity requires that the tensor acts on inner product spaces, and the current implementation requires InnerProductStyle(t) === EuclideanProduct().
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigh!(t). Note that the permuted tensor on which eigh! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
eig(t::AbstractTensor, (leftind, rightind)::Index2Tuple; kwargs...) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eig assumes that the linear map is not hermitian and returns type stable complex valued D and V tensors for both real and complex valued t. See eigh for hermitian linear maps
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eig!(t). Note that the permuted tensor on which eig! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
permute(t, (leftind, rightind)) * V = V * D
Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.
eigen(t::AbstractTensor, (leftind, rightind)::Index2Tuple; kwargs...) -> D, V
Compute eigenvalue factorization of tensor t as linear map from rightind to leftind.
If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigen!(t). Note that the permuted tensor on which eigen! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy
permute(t, (leftind, rightind)) * V = V * D
Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.
From categories to anyons: a travelogue
Kerstin Beer, Dmytro Bondarenko, Alexander Hahn, Maria Kalabakov, Nicole Knust, Laura Niermann, Tobias J. Osborne, Christin Schridde, Stefan Seckmeyer, Deniz E. Stiegemann, and Ramona Wolf
- [https://arxiv.org/abs/1811.06670](https://arxiv.org/abs/1811.06670)
Settings
This document was generated with Documenter.jl version 0.27.25 on Wednesday 1 November 2023. Using Julia version 1.9.3.