Symmetry sectors
Type hierarchy
The fundamental abstract supertype for symmetry sectors is Sector:
TensorKitSectors.Sector — Type
abstract type SectorAbstract 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:
unit(::Type{I}): unit element ofI. If there are multiple, implementallunits(::Type{I})instead.dual(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 timescappears ina ⊗ b, i.e. the multiplicityFusionStyle(::Type{I}):UniqueFusion(),SimpleFusion()orGenericFusion()BraidingStyle(::Type{I}):Bosonic(),Fermionic(),Anyonic(), ...Fsymbol(a::I, b::I, c::I, d::I, e::I, f::I): F-symbol: scalar (in case ofUniqueFusion/SimpleFusion) or matrix (in case ofGenericFusion)Rsymbol(a::I, b::I, c::I): R-symbol: scalar (in case ofUniqueFusion/SimpleFusion) or matrix (in case ofGenericFusion)
and optionally
dim(a::I): quantum dimension of sectorafrobenius_schur_indicator(a::I): Frobenius-Schur indicator ofa(1, 0, -1)frobenius_schur_phase(a::I): Frobenius-Schur phase ofa(±1)sectorscalartype(::Type{I}): scalar type of F- and R-symbolsBsymbol(a::I, b::I, c::I): B-symbol: scalar (in case ofUniqueFusion/SimpleFusion) or matrix (in case ofGenericFusion)twist(a::I)-> twist of sectora
Furthermore, iterate and Base.IteratorSize should be made to work for the singleton type SectorValues{I}.
To help with the implementation of ⊗(a::I, b::I) as an iterator, the provided struct type SectorProductIterator{I} can be used, which stores a and b and requires the implementation of Base.iterate(::SectorProductIterator{I}, state...).
Various concrete subtypes of Sector are provided within the TensorKitSectors library:
TensorKitSectors.Trivial — Type
struct Trivial <: Sector
Trivial()Singleton type to represent the trivial sector, i.e. the trivial representation of the trivial group. This is equivalent to Rep[ℤ₁], or the unit object of the category Vect of ordinary vector spaces.
TensorKitSectors.AbstractIrrep — Type
abstract type AbstractIrrep{G <: Group} <: SectorAbstract 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.
TensorKitSectors.ZNIrrep — Type
struct ZNIrrep{N} <: AbstractIrrep{ℤ{N}}
ZNIrrep{N}(n::Integer)
Irrep[ℤ{N}](n::Integer)Represents irreps of the group $ℤ_N$ for some value of N. For N equals 2, 3 or 4, ℤ{N} can be replaced by ℤ₂, ℤ₃, and ℤ₄. An arbitrary Integer n can be provided to the constructor, but only the value mod(n, N) is relevant.
The type of the stored integer (UInt8) requires N ≤ 128. Larger values of N should use the LargeZNIrrep instead. The constructor Irrep[ℤ{N}] should be preferred, as it will automatically select the most efficient storage type for a given value of N.
See also charge and modulus to extract the relevant data.
Fields
n::UInt8: the integer label of the irrep, moduloN.
TensorKitSectors.DNIrrep — Type
struct DNIrrep{N} <: AbstractIrrep{Dihedral{N}}
DNIrrep{N}(n::Integer, isodd::Bool=false)
Irrep[Dihedral{N}](n::Integer, isodd::Bool=false)Represents irreps of the dihedral group $D_N = Z_N ⋊ C$ ($Z_N$ and charge conjugation or reflection).
Properties
j::Int: the value of the $Z_N$ charge.isodd::Bool: the representation of charge conjugation.
Combined these take the values $+0, -0, 1, ..., (N - 1) / 2$ for odd $N$, and $+0, -0, 1, ..., N / 2 - 1, +(N/2), -(N/2)$ for even $N$, where the $+$ ($-$) refer to the even (odd) one-dimensional irreps, while the others are two-dimensional.
sourceTensorKitSectors.U1Irrep — Type
struct U1Irrep <: AbstractIrrep{U₁}
U1Irrep(charge::Real)
Irrep[U₁](charge::Real)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 of $U₁$ subgroups of $SU₂$, such as the $S^z$ of spin-1/2 system. Hence, the charge is stored as a HalfInt from the package HalfIntegers.jl, but can be entered as arbitrary Real. The sequence of the charges is: 0, 1/2, -1/2, 1, -1, ...
Fields
charge::HalfInt: the label of the irrep, which can be any half integer.
TensorKitSectors.SU2Irrep — Type
struct SU2Irrep <: AbstractIrrep{SU₂}
SU2Irrep(j::Real)
Irrep[SU₂](j::Real)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.
Fields
j::HalfInt: the label of the irrep, which can be any non-negative half integer.
TensorKitSectors.CU1Irrep — Type
struct CU1Irrep <: AbstractIrrep{CU₁}
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₂.
Fields
j::HalfInt: the value of the $U₁$ charge.s::Int: the representation of charge conjugation.
They can take values:
- if
j == 0,s = 0(trivial charge conjugation) ors = 1(non-trivial charge conjugation) - if
j > 0,s = 2(two-dimensional representation)
TensorKitSectors.AbstractGroupElement — Type
abstract type AbstractGroupElement{G <: Group} <: SectorAbstract supertype for sectors which corresponds to group elements of a group G.
Actual concrete implementations of those irreps can be obtained as Element[G], or via their actual name, which generically takes the form (asciiG)Element, i.e. the ASCII spelling of the group name followed by Element.
All group elements have FusionStyle equal to UniqueFusion(). Furthermore, the BraidingStyle is set to NoBraiding(), although this can be overridden by a concrete implementation of AbstractGroupElement.
For the fusion structure, a specific SomeGroupElement <: AbstractGroupElement{SomeGroup} should only implement the following methods
Base.:*(c1::GroupElement, c2::GroupElement) -> GroupElement
Base.one(::Type{GroupElement}) -> GroupElement
Base.inv(c::GroupElement) -> GroupElement
# and optionally
TensorKitSectors.cocycle(c1::GroupElement, c2::GroupElement, c3::GroupElement) -> NumberThe methods conj, dual, ⊗, Nsymbol, Fsymbol, dim, Asymbol, Bsymbol and frobenius_schur_phase will then be automatically defined. If no cocycle method is defined, the cocycle will be assumed to be trivial, i.e. equal to 1.
TensorKitSectors.ZNElement — Type
struct ZNElement{N, p} <: AbstractGroupElement{ℤ{N}}
ZNElement{N, p}(n::Integer)
GroupElement[ℤ{N}, p](n::Integer)Represents an element 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 ℤ₂, ℤ₃, ℤ₄. An arbitrary Integer n can be provided to the constructor, but only the value mod(n, N) is relevant. The second type parameter p should also be specified as an integer 0 <= p < N and specifies the 3-cocycle, which is then being given by
cocycle(a, b, c) = cispi(2 * p * a.n * (b.n + c.n - mod(b.n + c.n, N)) / N))If p is not specified, it defaults to 0, i.e. the trivial cocycle.
Fields
n::Int8: the integer label of the element, moduloN.
TensorKitSectors.FermionParity — Type
struct FermionParity <: Sector
FermionParity(isodd::Bool)Represents sectors with fermion parity. The fermion parity is a $ℤ₂$ quantum number that yields an additional sign when two odd fermions are exchanged, corresponding to a BraidingStyle that is Fermionic().
Fields
isodd::Bool: indicates whether the fermion parity is odd (true) or even (false).
See also: FermionNumber, FermionSpin
TensorKitSectors.FermionNumber — Type
const FermionNumber = U1Irrep ⊠ FermionParity
FermionNumber(a::Int)Represents the fermion number as the direct product of a $U₁$ irrep a and a fermion parity, with the restriction that the fermion parity is odd if and only if a is odd.
See also: U1Irrep, FermionParity
TensorKitSectors.FermionSpin — Type
const FermionSpin = SU2Irrep ⊠ FermionParity
FermionSpin(j::Real)Represents the fermion spin as the direct product of a $SU₂$ irrep j and a fermion parity, with the restriction that the fermion parity is odd if 2 * j is odd.
See also: SU2Irrep, FermionParity
TensorKitSectors.FibonacciAnyon — Type
struct FibonacciAnyon <: Sector
FibonacciAnyon(s::Symbol)Represents the anyons of the Fibonacci modular fusion category. It can take two values, corresponding to the trivial sector FibonacciAnyon(:I) and the non-trivial sector FibonacciAnyon(:τ) with fusion rules $τ ⊗ τ = 1 ⊕ τ$.
Fields
isunit::Bool: indicates whether the sector corresponds to the trivial anyon:I(true), or the non-trivial anyon:τ(false).
TensorKitSectors.IsingAnyon — Type
struct IsingAnyon <: Sector
IsingAnyon(s::Symbol)Represents the anyons of the Ising modular fusion category. It can take three values, corresponding to the trivial sector IsingAnyon(:I) and the non-trivial sectors IsingAnyon(:σ) and IsingAnyon(:ψ), with fusion rules $ψ ⊗ ψ = 1$, $σ ⊗ ψ = σ$, and $σ ⊗ σ = 1 ⊕ ψ$.
Fields
s::Symbol: the label of the represented anyon, which can be:I,:σ, or:ψ.
TensorKitSectors.PlanarTrivial — Type
struct PlanarTrivial <: Sector
PlanarTrivial()Represents a trivial anyon sector, i.e. a trivial sector without braiding. This is mostly useful for testing.
sourceTensorKitSectors.IsingBimodule — Type
struct IsingBimodule <: SectorType to represent the simple objects in the Ising category reinterpreted as a bimodule category composed of two copies of the category 𝒞 = 𝒟 = Irrep[ℤ₂], the two simple objects of which can be identified with the Ising anyons {I, ψ}, and the bimodule categories ℳ = ℳᵒᵖ = Vec, with a single simple object that can be identified with the Ising anyon σ. This constitutes the easiest example of a multifusion category and is implemented here for testing purposes and to illustrate how to implement such categories in TensorKitSectors.jl.
sourceTensorKitSectors.TimeReversed — Type
struct TimeReversed{I <: Sector}
TimeReversed(a::I) where {I <: Sector}Represents the time-reversed version of the sector I, i.e. the sector with the same fusion rules and F-symbols, but with the inverse braiding.
TensorKitSectors.ProductSector — Type
struct ProductSector{T <: SectorTuple}
ProductSector((s₁, s₂, ...))Represents the Deligne tensor product of sectors. The type parameter T is a tuple of the component sectors. The recommended way to construct a ProductSector is using the deligneproduct (⊠) operator on the components.
Several more concrete sector types can be found in other packages such as SUNRepresentations.jl, CategoryData.jl, QWignerSymbols.jl, ...:
Some of these types are parameterized by a type parameter that represents a group. We therefore also provide a number of types to represent groups:
TensorKitSectors.Group — Type
abstract type GroupAbstract supertype for representing different types of groups. Groups can be used to define Sector subtypes, either via their irreducible representations, or via their group elements, and typically appear as type parameter. As such, they are not meant to be instantiated and are defined as abstract types.
sourceTensorKitSectors.AbelianGroup — Type
abstract type AbelianGroup <: GroupAbstract supertype for representing different types of Abelian groups. Abelian groups have both irreps and group elements that have several simplified properties, that can be defined in general.
sourceTensorKitSectors.Cyclic — Type
abstract type Cyclic{N} <: AbelianGroupType to represent the cyclic group of order N, i.e. the multiplicative group of roots of unity of order N, which is a discrete abelian group. The cyclic group of order N is isomorphic to the additive group ℤ{N}, and we define the latter as a type alias const ℤ{N} = Cyclic{N}.
TensorKitSectors.U₁ — Type
abstract type U₁ <: AbelianGroupType to represent the group $U(1)$ of complex numbers of unit modulus, which is a compact Abelian Lie group.
sourceTensorKitSectors.CU₁ — Type
abstract type CU₁ <: GroupType to represent the group of U₁ in combination with charge conjugation, i.e. the group generated by U₁ and an additional element that acts as complex conjugation on U₁. This group is isomorphic to the orthogonal group O₂ of real orthogonal 2×2 matrices, and can be seen as the semidirect product U₁ ⋊ ℤ₂. This is a compact non-Abelian group.
sourceTensorKitSectors.SU — Type
abstract type SU{N} <: GroupType to represent the special unitary group $SU(N)$, which is a compact non-Abelian Lie group.
sourceTensorKitSectors.Dihedral — Type
abstract type Dihedral{N} <: GroupType to represent the dihedral group of order 2N, which is the symmetry group of a regular polygon with N sides, and is a discrete non-Abelian group.
The following types are used to characterise different properties of the different types of sectors:
TensorKitSectors.FusionStyle — Type
abstract type FusionStyle
FusionStyle(::Sector)
FusionStyle(I::Type{<:Sector})Trait to describe the 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}.
TensorKitSectors.BraidingStyle — Type
abstract type BradingStyle
BraidingStyle(::Sector) -> ::BraidingStyle
BraidingStyle(I::Type{<:Sector}) -> ::BraidingStyleReturn the type of braiding and twist behavior of sectors of type I, which can be either
NoBraiding(): no braiding structureBosonic(): symmetric braiding with trivial twist (i.e. identity)Fermionic(): symmetric braiding with non-trivial twist (squares to identity)Anyonic(): general $R^{ab}_c$ phase or matrix (depending onSimpleFusionorGenericFusionfusion) 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.
TensorKitSectors.UnitStyle — Type
abstract type UnitStyle
UnitStyle(::Sector)
UnitStyle(I::Type{<:Sector})Trait to describe the semisimplicity of the unit sector of type I. This can be either
SimpleUnit(): the unit is simple (e.g. fusion categories);GenericUnit(): the unit is semisimple.
Finally, the following auxiliary types are defined to facilitate the implementation of some of the methods on sectors:
TensorKitSectors.SectorValues — Type
struct SectorValues{I <: Sector}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 valuesBase.IteratorSize(::Type{SectorValues{I}}):HasLength(),SizeUnknown()orIsInfinite()depending on whether the number of values of typeIis finite (and sufficiently small) or infinite; for a large number of values,SizeUnknown()is recommended because this will trigger the use ofGenericGradedSpace.
If IteratorSize(I) == HasLength(), also the following must be implemented:
Base.length(::SectorValues{I}): the number of different valuesBase.getindex(::SectorValues{I}, i::Int): a mapping between an indexiand an instance ofI. A fallback implementation exists that returns theith value of theSectorValuesiterator.findindex(::SectorValues{I}, c::I): reverse mapping between a valuec::Iand an indexi::Integer ∈ 1:length(values(I)). A fallback implementation exists that linearly searches through theSectorValuesiterator.
TensorKitSectors.SectorProductIterator — Type
struct SectorProductIterator{I <: Sector}
SectorProductIterator(a::I, b::I) where {I <: Sector}Custom iterator to represent the (unique) fusion outputs of $a ⊗ b$.
Custom sectors that aim to use this have to provide the following functionality:
Base.iterate(::SectorProductIterator{I}, state...) where {I <: Sector}: iterate over the fusion outputs ofa ⊗ b
If desired and it is possible to easily compute the number of unique fusion outputs, it is also possible to define Base.IteratorSize(::Type{SectorProductIterator{I}}) = Base.HasLength(), in which case Base.length(::SectorProductIterator{I}) has to be implemented.
See also ⊗.
Useful constants
The following constants are defined to facilitate obtaining the type associated with the group elements or the irreducible representations of a given group:
TensorKitSectors.Irrep — Constant
const IrrepA 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.
TensorKitSectors.GroupElement — Constant
const GroupElementA constant of a singleton type used as GroupElement[G] or GroupElement[G, ω] with G <: Group a type of group, to construct or obtain a concrete subtype of AbstractElement{G} that implements the data structure used to represent elements of the group G, possibly with a second argument ω that specifies the associated 3-cocycle.
Methods for characterizing and manipulating Sector objects
The following methods can be used to obtain properties such as topological data of sector objects, or to manipulate them or create related sectors:
TensorKitSectors.unit — Function
unit(::Sector) -> Sector
unit(::Type{<:Sector}) -> SectorReturn the unit element of this type of sector, provided it is unique.
sourceTensorKitSectors.isunit — Function
TensorKitSectors.leftunit — Function
TensorKitSectors.rightunit — Function
TensorKitSectors.allunits — Function
allunits(I::Type{<:Sector}) -> Tuple{I}Return a tuple with all units of the sector type I. For fusion categories, this will contain only one element.
TensorKitSectors.dual — Method
dual(a::Sector) -> SectorReturn the dual label of a, i.e. the unique label ā = dual(a) such that Nsymbol(a, ā, leftunit(a)) == 1 and Nsymbol(ā, a, rightunit(a)) == 1.
TensorKitSectors.Nsymbol — Function
Nsymbol(a::I, b::I, c::I) where {I <: Sector} -> IntegerReturn 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().
TensorKitSectors.:⊗ — Function
⊗(a::I, b::I...) where {I <: Sector}
otimes(a::I, b::I...) where {I <: Sector}Return an iterable of elements of c::I that appear in the fusion product a ⊗ b.
Note that every element c should appear at most once, fusion degeneracies (if FusionStyle(I) == GenericFusion()) should be accessed via Nsymbol(a, b, c).
TensorKitSectors.Fsymbol — Function
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-<-κ
∨
cIf 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)).
TensorKitSectors.Rsymbol — Function
Rsymbol(a::I, b::I, c::I) where {I <: Sector}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 aIf 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).
TensorKitSectors.Bsymbol — Function
Bsymbol(a::I, b::I, c::I) where {I <: Sector}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).
TensorKitSectors.dim — Method
TensorKitSectors.frobenius_schur_phase — Function
frobenius_schur_phase(a::Sector)Return the Frobenius-Schur phase $κₐ$ of a sector $a$, which is a complex phase that appears in the context of bending lines and is obtained from $F^{a a̅ a}_a$. When a == dual(a), it is restricted to $κₐ ∈ \{1, -1\}$ and coincides with the group-theoretic version frobenius_schur_indicator. When a != dual(a), the value of $κₐ$ can be gauged to be 1, though is not required to be.
TensorKitSectors.frobenius_schur_indicator — Function
frobenius_schur_indicator(a::Sector)Return the Frobenius-Schur indicator of a sector $νₐ ∈ \{1, 0, -1\}$, which distinguishes between real, complex and quaternionic representations.
See also frobenius_schur_phase for the category-theoretic version that appears in the context of line bending.
TensorKitSectors.twist — Method
Base.isreal — Method
isreal(::Type{<:Sector}) -> BoolReturn whether the topological data (Fsymbol, Rsymbol) of the sector is real or not (in which case it is complex).
sourceTensorKitSectors.sectorscalartype — Function
sectorscalartype(I::Type{<:Sector}) -> TypeReturn the scalar type of the topological data (Fsymbol, Rsymbol) of the sector I.
TensorKitSectors.deligneproduct — Method
⊠(s₁::Sector, s₂::Sector)
deligneproduct(s₁::Sector, s₂::Sector)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 \boxtimes+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₂].
We have also the following methods that are specific to certain types of sectors and serve as accessors to their fields:
TensorKitSectors.charge — Function
TensorKitSectors.modulus — Function
modulus(c::ZNIrrep{N}) -> N
modulus(::Type{<:ZNIrrep{N}}) -> NThe order of the cyclic group, or the modulus of the charge labels.
sourceFurthermore, we also have one specific method acting on groups, represented as types
TensorKitSectors.:× — Function
×(G::Vararg{Type{<:Group}}) -> ProductGroup{Tuple{G...}}
times(G::Vararg{Type{<:Group}}) -> ProductGroup{Tuple{G...}}Construct the direct product of a (list of) groups.
sourceBecause we sometimes want to customize the string representation of our sector types, we also have the following method:
TensorKitSectors.type_repr — Function
type_repr(T::Type)Return a string representation of the type T, which is used to modify the default way in which Sector subtypes are displayed in other objects that depend on them.
Finally, we provide functionality to compile all revelant methods for a sector: