Grassmann.jl Library

AbstractTensors.BivectorType
Bivector{V,T} <: TensorGraded{V,2,T}

Graded bivector elements of a Manifold instance V with scalar field T.

AbstractTensors.GradedVectorType
GradedVector{V,T} <: TensorGraded{V,1,T}

Graded vector elements of a Manifold instance V with scalar field T.

AbstractTensors.ManifoldType
Manifold{V,T} <: TensorAlgebra{V,T}

Basis parameter locally homeomorphic to V::Submanifold{M} T-module product topology.

AbstractTensors.ScalarType
Scalar{V,T} <: TensorGraded{V,0,T}

Graded scalar elements of a Manifold instance V with scalar field T.

AbstractTensors.TensorGradedType
TensorGraded{V,G,T} <: Manifold{V,T} <: TensorAlgebra

Grade G elements of a Manifold instance V with scalar field T.

AbstractTensors.TensorMixedType
TensorMixed{V,T} <: TensorAlgebra{V,T}

Elements of Manifold instance V having non-homogenous grade with scalar field T.

AbstractTensors.TensorTermType
TensorTerm{V,G,T} <: TensorGraded{V,G,T}

Single coefficient for grade G of a Manifold instance V with scalar field T.

AbstractTensors.TrivectorType
Trivector{V,T} <: TensorGraded{V,3,T}

Graded trivector elements of a Manifold instance V with scalar field T.

AbstractTensors.antisandwichMethod
antisandwich(x::TensorAlgebra,R::TensorAlgebra)

Defined as complementleft(complementright(R)>>>complementright(x)).

AbstractTensors.bivectorFunction
bivector(::TensorAlgebra)

Return the bivector (rank 2) part of any TensorAlgebra element.

AbstractTensors.coabsFunction
coabs(t::TensorAlgebra)

Complemented abs defined as complementleft(abs(complementright(t))).

AbstractTensors.coabs2Function
coabs2(t::TensorAlgebra)

Complemented abs2 defined as complementleft(abs2(complementright(t))).

AbstractTensors.cocbrtFunction
cocbrt(t::TensorAlgebra)

Complemented cbrt defined as complementleft(cbrt(complementright(t))).

AbstractTensors.cocosFunction
cocos(t::TensorAlgebra)

Complemented cos defined as complementleft(cos(complementright(t))).

AbstractTensors.cocoshFunction
cocosh(t::TensorAlgebra)

Complemented cosh defined as complementleft(cosh(complementright(t))).

AbstractTensors.coexpFunction
coexp(t::TensorAlgebra)

Complemented exp defined as complementleft(exp(complementright(t))).

AbstractTensors.coinvFunction
coinv(t::TensorAlgebra)

Complemented inv defined as complementleft(inv(complementright(t))).

AbstractTensors.cologFunction
colog(t::TensorAlgebra)

Complemented log defined as complementleft(log(complementright(t))).

AbstractTensors.cosandwichMethod
cosandwich(x::TensorAlgebra,R::TensorAlgebra)

Defined as complementleft(sandwich(complementright(x),complementright(R))).

AbstractTensors.cosinFunction
cosin(t::TensorAlgebra)

Complemented sin defined as complementleft(sin(complementright(t))).

AbstractTensors.cosinhFunction
cosinh(t::TensorAlgebra)

Complemented sinh defined as complementleft(sinh(complementright(t))).

AbstractTensors.cosqrtFunction
cosqrt(t::TensorAlgebra)

Complemented sqrt defined as complementleft(sqrt(complementright(t))).

AbstractTensors.cotanFunction
cotan(t::TensorAlgebra)

Complemented tan defined as complementleft(tan(complementright(t))).

AbstractTensors.cotanhFunction
cotanh(t::TensorAlgebra)

Complemented tanh defined as complementleft(tanh(complementright(t))).

AbstractTensors.counitMethod
counit(t::TensorAlgebra)

Pseudo-normalization defined as unitize(t) = t/value(coabs(t)).

AbstractTensors.gdimsMethod
gdims(t::TensorGraded{V,G})

Dimensionality of the grade G of V for that TensorAlgebra.

AbstractTensors.geomabsMethod
geomabs(t::TensorAlgebra)

Geometric norm defined as geomabs(t) = abs(t) + coabs(t).

AbstractTensors.mdimsMethod
mdims(t::TensorAlgebra{V})

Dimensionality of the pseudoscalar V of that TensorAlgebra.

AbstractTensors.pseudoabsFunction
pseudoabs(t::TensorAlgebra)

Complemented abs defined as complementleft(abs(complementright(t))).

AbstractTensors.pseudoabs2Function
pseudoabs2(t::TensorAlgebra)

Complemented abs2 defined as complementleft(abs2(complementright(t))).

AbstractTensors.pseudocbrtFunction
pseudocbrt(t::TensorAlgebra)

Complemented cbrt defined as complementleft(cbrt(complementright(t))).

AbstractTensors.pseudocosFunction
pseudocos(t::TensorAlgebra)

Complemented cos defined as complementleft(cos(complementright(t))).

AbstractTensors.pseudocoshFunction
pseudocosh(t::TensorAlgebra)

Complemented cosh defined as complementleft(cosh(complementright(t))).

AbstractTensors.pseudoexpFunction
pseudoexp(t::TensorAlgebra)

Complemented exp defined as complementleft(exp(complementright(t))).

AbstractTensors.pseudoinvFunction
pseudoinv(t::TensorAlgebra)

Complemented inv defined as complementleft(inv(complementright(t))).

AbstractTensors.pseudologFunction
pseudolog(t::TensorAlgebra)

Complemented log defined as complementleft(log(complementright(t))).

AbstractTensors.pseudoscalarMethod
pseudoscalar(::TensorAlgebra), volume(::TensorAlgebra)

Return the pseudoscalar (full rank) part of any TensorAlgebra element.

AbstractTensors.pseudosinFunction
pseudosin(t::TensorAlgebra)

Complemented sin defined as complementleft(sin(complementright(t))).

AbstractTensors.pseudosinhFunction
pseudosinh(t::TensorAlgebra)

Complemented sinh defined as complementleft(sinh(complementright(t))).

AbstractTensors.pseudosqrtFunction
pseudosqrt(t::TensorAlgebra)

Complemented sqrt defined as complementleft(sqrt(complementright(t))).

AbstractTensors.pseudotanFunction
pseudotan(t::TensorAlgebra)

Complemented tan defined as complementleft(tan(complementright(t))).

AbstractTensors.pseudotanhFunction
pseudotanh(t::TensorAlgebra)

Complemented tanh defined as complementleft(tanh(complementright(t))).

AbstractTensors.scalarFunction
scalar(::TensorAlgebra)

Return the scalar (rank 0) part of any TensorAlgebra element.

AbstractTensors.tdimsMethod
tdims(t::TensorAlgebra{V})

Dimensionality of the superalgebra of V for that TensorAlgebra.

AbstractTensors.trivectorFunction
trivector(::TensorAlgebra)

Return the trivector (rank 3) part of any TensorAlgebra element.

AbstractTensors.unitnormMethod
unitnorm(t::TensorAlgebra)

Geometric normalization defined as unitnorm(t) = t/norm(geomabs(t)).

AbstractTensors.valueMethod
value(::TensorAlgebra)

Returns the internal Values representation of a TensorAlgebra element.

AbstractTensors.valuetypeMethod
valuetype(t::TensorAlgebra{V,T}) where {V,T} = T

Returns type of a TensorAlgebra element value's internal representation.

AbstractTensors.vectorFunction
vector(::TensorAlgebra)

Return the vector (rank 1) part of any TensorAlgebra element.

LinearAlgebra.rankMethod
rank(::Manifold{n})

Dimensionality n of the Manifold subspace representation.

AbstractTensors.@coMacro
@co fun(args...)

Use the macro @co to make a pseudoscalar complement variant of any functions:

julia> @co myfun(x)
comyfun (generic function with 1 method)

Now comyfun(x) = complementleft(myfun(complementright(x))) is defined.

julia> @co myproduct(a,b)
comyproduct (generic function with 1 method)

Now comyproduct(a,b) = complementleft(myproduct(!a,!b)) is defined.

AbstractTensors.@pseudoMacro
@pseudo fun(args...)

Use the macro @pseudo to make a pseudoscalar complement variant of any functions:

julia> @pseudo myfun(x)
pseudomyfun (generic function with 1 method)

Now pseudomyfun(x) = complementleft(myfun(complementright(x))) is defined.

julia> @pseudo myproduct(a,b)
pseudomyproduct (generic function with 1 method)

Now pseudomyproduct(a,b) = complementleft(myproduct(!a,!b)) is defined.

DirectSum.BasisType
DirectSum.Basis{V} <: SubAlgebra{V} <: TensorAlgebra{V}

Grassmann basis container with cache of Submanifold elements and their Symbol names.

DirectSum.ExtendedBasisType
DirectSum.ExtendedBasis{V} <: SubAlgebra{V} <: TensorAlgebra{V}

Grassmann basis container without a dedicated Submanifold cache (only lazy caching).

DirectSum.InfinityType
Infinity{V} <: TensorGraded{V,0} <: TensorAlgebra{V}

Infinite quantity Infinity of the Grassmann algebra over V.

DirectSum.OneType
One{V} <: TensorGraded{V,0} <: TensorAlgebra{V}

Unit quantity One of the Grassmann algebra over V.

DirectSum.SingleType
Single{V,G,B,T} <: TensorTerm{V,G,T} <: TensorGraded{V,G,T}

Single type with pseudoscalar V::Manifold, grade/rank G::Int, B::Submanifold{V,G}, field T::Type.

DirectSum.SparseBasisType
DirectSum.SparseBasis{V} <: SubAlgebra{V} <: TensorAlgebra{V}

Grassmann basis with sparse cache of Submanifold{G,V} elements and their Symbol names.

DirectSum.SubmanifoldType
Submanifold{V,G,B} <: TensorGraded{V,G} <: Manifold{G}

Basis type with pseudoscalar V::Manifold, grade/rank G::Int, bits B::UInt64.

DirectSum.TensorBundleType
TensorBundle{n,ℙ,g,ν,μ} <: Manifold{n}

A manifold representing the space of tensors with a given metric and basis for projection.

Let n be the rank of a Manifold{n}. The type TensorBundle{n,ℙ,g,ν,μ} uses byte-encoded data available at pre-compilation, where specifies the basis for up and down projection, g is a bilinear form that specifies the metric of the space, and μ is an integer specifying the order of the tangent bundle (i.e. multiplicity limit of Leibniz-Taylor monomials). Lastly, ν is the number of tangent variables.

DirectSum.ZeroType
Zero{V} <: TensorGraded{V,0} <: TensorAlgebra{V}

Null quantity Zero of the Grassmann algebra over V.

AbstractTensors.cliffordFunction
clifford(ω::TensorAlgebra)

Clifford conjugate of an element: clifford(ω) = involute(reverse(ω))

Base.:~Method
reverse(ω::TensorAlgebra)

Reverse of an element: ~ω = (-1)^(grade(ω)(grade(ω)-1)/2)ω

Base.conjFunction
~(ω::TensorAlgebra)

Reverse of an element: ~ω = (-1)^(grade(ω)(grade(ω)-1)/2)ω

Base.imagMethod
imag(ω::TensorAlgebra)

The imag part (ω-(~ω))/2 is defined by abs2(imag(ω)) == -(imag(ω)^2).

Base.realMethod

real(ω::TensorAlgebra)

The real part (ω+(~ω))/2 is defined by abs2(real(ω)) == real(ω)^2.

DirectSum.allocFunction
alloc(V::Manifold,:V,"v","w","∂","ϵ")

Generates Basis declaration having Manifold specified by V. The first argument provides pseudoscalar specifications, the second argument is the variable name for the Manifold, and the third and fourth argument are variable prefixes of the Submanifold vector names (and covector basis names).

DirectSum.getbasisMethod
getbasis(V::Manifold,v)

Fetch a specific Submanifold{G,V} element from an optimal SubAlgebra{V} selection.

DirectSum.nameindexMethod
nameindex(V::Int) -> Int

Returns a default name index for an integer input.

DirectSum.nameindexMethod
nameindex(a::NTuple{4, String}) -> Int

Get a unique index for a set of name options, based on a cached result.

DirectSum.nameindexMethod
nameindex(V::T) -> Int

Returns the name index for V when V is a TensorBundle or Manifold.

DirectSum.namelistMethod
namelist(V) -> String

Get the cached name list for the options specified in V.

DirectSum.pseudocliffordFunction
pseudoclifford(ω::TensorAlgebra)

Pseudo-clifford conjugate element: pseudoclifford(ω) = pseudoinvolute(pseudoreverse(ω))

DirectSum.pseudoinvoluteFunction
pseudoinvolute(ω::TensorAlgebra)

Anti-involute of an element: ~ω = (-1)^pseudograde(ω)*ω

DirectSum.pseudoreverseFunction
pseudoreverse(ω::TensorAlgebra)

Anti-reverse of an element: ~ω = (-1)^(pseudograde(ω)(pseudograde(ω)-1)/2)ω

DirectSum.tensorhashFunction
tensorhash(d, o; c=0, C=0) -> Int

Compute a unique integer representation for a set of options.

Arguments:

  • d: Number of tangent variables (corresponds to ν).
  • o: Order of Leibniz-Taylor monomials multiplicity limit (corresponds to μ).
  • c: Coefficient used to specify the metric (corresponds to g). Defaults to 0.
  • C: Additional coefficient. Defaults to 0.

Returns

A unique integer representation that can be used as a cache key.

DirectSum.@basisMacro
@basis

Generates Submanifold elements having Manifold specified by V. As a result of this macro, all of the Submanifold{V,G} elements generated by that TensorBundle become available in the local workspace with the specified naming. The first argument provides pseudoscalar specifications, the second argument is the variable name for the Manifold, and the third and fourth argument are variable prefixes of the Submanifold vector names (and covector basis names). Default for @basis M is @basis M V v w ∂ ϵ.

DirectSum.@dualbasisMacro
@dualbasis

Generates Submanifold elements having Manifold specified by V'. As a result of this macro, all of the Submanifold{V',G} elements generated by that TensorBundle become available in the local workspace with the specified naming. The first argument provides pseudoscalar specifications, the second argument is the variable name for the dual Manifold, and the third and fourth argument are variable prefixes of the Submanifold covector names (and tensor field basis names). Default for @dualbasis M is @dualbasis M VV w ϵ.

DirectSum.@mixedbasisMacro
@mixedbasis

Generates Submanifold elements having Manifold specified by V⊕V'. As a result of this macro, all of the Submanifold{V⊕V',G} elements generated by that TensorBundle become available in the local workspace with the specified naming. The first argument provides pseudoscalar specifications, the second argument is the variable name for the Manifold, and the third and fourth argument are variable prefixes of the Submanifold vector names (and covector basis names). Default for @mixedbasis M is @mixedbasis M V v w ∂ ϵ.

Grassmann.AbstractSpinorType
AbstractSpinor{V,T} <: TensorMixed{V,T} <: TensorAlgebra{V,T}

Elements of TensorAlgebra having non-homogenous grade being a spinor in the abstract.

Grassmann.ChainMethod
Chain{V,G,T} <: TensorGraded{V,G,T} <: TensorAlgebra{V,T}

Chain type with pseudoscalar V::Manifold, grade/rank G::Int, scalar field T::Type.

Grassmann.ChainBundleType
ChainBundle{V,G,T,P} <: Manifold{V,T} <: TensorAlgebra{V,T}

Subsets of a bundle cross-section over a Manifold topology.

Grassmann.CoSpinorMethod
CoSpinor{V,T} <: AbstractSpinor{V,T} <: TensorAlgebra{V,T}

PsuedoSpinor (odd grade) type with pseudoscalar V::Manifold and scalar T::Type.

Grassmann.CoupleType
Couple{V,B,T} <: AbstractSpinor{V,T} <: TensorAlgebra{V,T}

Pair of values with V::Manifold, basis B::Submanifold, scalar field of T::Type.

Grassmann.MultivectorMethod
Multivector{V,T} <: TensorMixed{V,T} <: TensorAlgebra{V,T}

Chain type with pseudoscalar V::Manifold and scalar field T::Type.

Grassmann.PhasorType
Phasor{V,B,T} <: AbstractSpinor{V,T} <: TensorAlgebra{V,T}

Magnitude/phase angle with V::Manifold, frequency B::Type, and amplitude T::Type.

Grassmann.PseudoCoupleType
PseudoCouple{V,B,T} <: AbstractSpinor{V,T} <: TensorAlgebra{V,T}

Pair of values with V::Manifold, basis B::Submanifold, pseudoscalar of T::Type.

Grassmann.SpinorMethod
Spinor{V,T} <: AbstractSpinor{V,T} <: TensorAlgebra{V,T}

Spinor (even grade) type with pseudoscalar V::Manifold and scalar field T::Type.

AbstractLattices.:∧Method
∧(ω::TensorAlgebra,η::TensorAlgebra)

Exterior product as defined by the anti-symmetric quotient Λ≡⊗/~

AbstractLattices.:∨Method
∨(ω::TensorAlgebra,η::TensorAlgebra)

Regressive product as defined by the DeMorgan's law: ∨(ω...) = ⋆⁻¹(∧(⋆.(ω)...))

AbstractTensors.:∗Function
∗(ω::TensorAlgebra,η::TensorAlgebra)

Reversed geometric product: ω∗η = (~ω)*η

AbstractTensors.:⊘Function
⊘(ω::TensorAlgebra,η::TensorAlgebra)

General sandwich product: ω⊘η = reverse(η)⊖ω⊖involute(η)

For normalized even grade η it is ω⊘η = (~η)⊖ω⊖η

AbstractTensors.:⊙Method
⊙(ω::TensorAlgebra,η::TensorAlgebra)

Symmetrization projection: ⊙(ω...) = ∑(∏(σ.(ω)...))/factorial(length(ω))

AbstractTensors.:⊠Method
⊠(ω::TensorAlgebra,η::TensorAlgebra)

Anti-symmetrization projection: ⊠(ω...) = ∑(∏(πσ.(ω)...))/factorial(length(ω))

AbstractTensors.:⟑Method
*(ω::TensorAlgebra,η::TensorAlgebra)

Geometric algebraic product: ω⊖η = (-1)ᵖdet(ω∩η)⊗(Λ(ω⊖η)∪L(ω⊕η))

AbstractTensors.contractionMethod
contraction(ω::TensorAlgebra,η::TensorAlgebra)

Interior (right) contraction product: ω⋅η = ω∨⋆η

AbstractTensors.evenMethod
even(t)

Selects the even part (t+involute(t))/2 and is defined by even grade.

AbstractTensors.oddMethod
odd(t)

Selects the odd part (t-involute(t))/2 and is defined by odd grade.

Base.:&Function
∨(ω::TensorAlgebra,η::TensorAlgebra)

Regressive product as defined by the DeMorgan's law: ∨(ω...) = ⋆⁻¹(∧(⋆.(ω)...))

Base.:>>>Function
>>>(ω::TensorAlgebra,η::TensorAlgebra)

Traditional sandwich product: ω>>>η = ω⊖η⊖clifford(ω)

For normalized even grade η it is ω>>>η = ω⊖η⊖(~ω)

Grassmann.:↑Function
↑(ω::TensorAlgebra{V}) where V # project

Canonical up-project operation from the space V, based on either Euclidean projective geometry, or the Riemann sphere, or conformal geometric algebra, or potentially other future canonical specifications. Optional arguments expose lower-level building blocks: ↑(ω,b) — use an explicit projective basis element b, or in the CGA context ↑(ω,p,m) — parameterise the conformal split with point-like part p and Minkowski part m. See also for the inverse down-reject operation.

Grassmann.:↓Function
↓(ω::TensorAlgebra{V}) where V # reject

Canonical down-reject operation from the space V, based on either Euclidean projective geometry, or the Riemann sphere, or conformal geometric algebra, or potentially other future canonical specifications. Optional arguments expose lower-level building blocks: ↓(ω,b) — use an explicit projective basis element b, or in the CGA context ↓(ω,p,m) — parameterise the conformal split with point-like part p and Minkowski part m. See also for the inverse up-project operation.

Grassmann.bettiMethod
betti(::TensorAlgebra)

Compute combinatoric Betti numbers based on the count_gdims and boundary_rank methods.

Grassmann.boundary_nullMethod
boundary_null(t)

Dimension of the nullspace (kernel) of the combinatorial boundary operator on t, which can be an element of the TensorAlgebra. For each grade k the function lists Values with entries dₖ₊₁ - rₖ, where the computation follows from d = count_gdims(t) and r = boundary_rank(t) utility methods. Intuitively, this measures how many k‑chains are boundaries of (k+1)‑chains and therefore vanish under the operation.

Grassmann.boundary_rankFunction
boundary_rank(t)

Returns the rank of the combinatorial boundary operator evaluated on t, which can be an element of the TensorAlgebra. Internally the function counts how many graded dimensions of t remain after applying the boundary operator on it, minimized to not exceed the original graded space. This results in a list of Values whose entries rₖ satisfy 0 ≤ rₖ ≤ dₖ computed with d = count_gdims(t) as utility methods.

Grassmann.cayleyFunction
cayley(V,op=*)

Compute the cayley table with op(a,b) for each Submanifold basis of V.

Grassmann.compactFunction
compactio(::Bool)

Toggles compact numbers for extended Grassmann elements (enabled by default)

Grassmann.projectFunction
↑(ω::TensorAlgebra{V}) where V # project

Canonical up-project operation from the space V, based on either Euclidean projective geometry, or the Riemann sphere, or conformal geometric algebra, or potentially other future canonical specifications. Optional arguments expose lower-level building blocks: ↑(ω,b) — use an explicit projective basis element b, or in the CGA context ↑(ω,p,m) — parameterise the conformal split with point-like part p and Minkowski part m. See also for the inverse down-reject operation.

Grassmann.rejectFunction
↓(ω::TensorAlgebra{V}) where V # reject

Canonical down-reject operation from the space V, based on either Euclidean projective geometry, or the Riemann sphere, or conformal geometric algebra, or potentially other future canonical specifications. Optional arguments expose lower-level building blocks: ↓(ω,b) — use an explicit projective basis element b, or in the CGA context ↓(ω,p,m) — parameterise the conformal split with point-like part p and Minkowski part m. See also for the inverse up-project operation.

LinearAlgebra.crossFunction
cross(ω::TensorAlgebra,η::TensorAlgebra)

Cross product: ω×η = ⋆(ω∧η)

LinearAlgebra.dotFunction
dot(ω::TensorAlgebra,η::TensorAlgebra)

Interior (right) contraction product: ω⋅η = ω∨⋆η

Leibniz.laplacianConstant
Δ = ∇^2 = ∂ₖ²v # laplacian::Laplacian

Abstract laplacian as second-order Derivation{Bool,2} is Laplacian operator dispatch.

Leibniz.nablaConstant
∇ = ∂ₖvₖ # nabla::Nabla

Abstract nabla as first-order Derivation{Bool,1} is Nabla operator dispatch.

Leibniz.ΔConstant
Δ = ∇^2 = ∂ₖ²v # laplacian::Laplacian

Abstract laplacian as second-order Derivation{Bool,2} is Laplacian operator dispatch.

Leibniz.∇Constant
∇ = ∂ₖvₖ # nabla::Nabla

Abstract nabla as first-order Derivation{Bool,1} is Nabla operator dispatch.

Leibniz.LaplacianType
Δ = ∇^2 = ∂ₖ²v # laplacian::Laplacian

Abstract laplacian as second-order Derivation{Bool,2} is Laplacian operator dispatch.

Leibniz.NablaType
∇ = ∂ₖvₖ # nabla::Nabla

Abstract nabla as first-order Derivation{Bool,1} is Nabla operator dispatch.

Leibniz.bit2intMethod
bit2int(b::BitVector) -> UInt

Takes a BitVector and converts it to a corresponding unsigned integer representation.

Leibniz.digitsfastMethod
digitsfast(b,N::Int) -> Values{N+1,Int}

Calculates the digits of a binary number b up to length N.

Leibniz.gdimsallMethod
gdimsall(N) -> Values{N+1,Int}

Returns Values{N+1,Int} representing binomial coefficients from 0 to N.

Leibniz.gdimsevenMethod
gdimseven(N) -> Values

Returns Values representing binomial coefficients for even values from 0 to N.

Leibniz.gdimsoddMethod
gdimsodd(N) -> Values

Returns Values representing binomial coefficients for odd values from 1 to N.

Leibniz.getbasisMethod
getbasis(V::Manifold,v)

Fetch a specific SubManifold{G,V} element from an optimal SubAlgebra{V} selection.

Leibniz.indexbitsMethod
indexbits(N,indices) -> BitVector

Converts a list of indices into a BitVector of length N, used to create an array where only the specified indices are set to true.

Leibniz.indicesMethod
indices(b::UInt) -> Vector
indices(b::UInt,N::UInt) -> Vector

Computes the indices at which a binary number b has bits equal to 1. The N argument (optional) specifies the length of the binary representation.

Leibniz.intlogMethod
intlog(M::Integer) -> Int

This function computes the integer 2-logarithm of an integer.

Leibniz.χMethod
χ(::TensorAlgebra)

Compute the Euler characteristic χ = ∑ₚ(-1)ᵖbₚ.