API of the TensND library
TensND.Basis
— TypeBasis(v::AbstractMatrix{T}, ::Val{:cov})
Basis{dim, T<:Number}()
Basis(θ::T<:Number, ϕ::T<:Number, ψ::T<:Number)
Basis built from a square matrix v
where columns correspond either to
- primal vectors ie
eᵢ=v[:,i]
ifvar=:cov
as by default - dual vectors ie
eⁱ=v[:,i]
ifvar=:cont
.
Basis without any argument refers to the canonical basis (CanonicalBasis
) in Rᵈⁱᵐ
(by default dim=3
and T=Sym
)
Basis can also be built from Euler angles (RotatedBasis
) θ
in 2D and (θ, ϕ, ψ)
in 3D
The attributes of this object can be obtained by
vecbasis(ℬ, :cov)
: square matrix defining the primal basiseᵢ=e[:,i]
vecbasis(ℬ, :cont)
: square matrix defining the dual basiseⁱ=E[:,i]
metric(ℬ, :cov)
: square matrix defining the covariant components of the metric tensorgᵢⱼ=eᵢ⋅eⱼ=g[i,j]
metric(ℬ, :cont)
: square matrix defining the contravariant components of the metric tensorgⁱʲ=eⁱ⋅eʲ=gⁱʲ[i,j]
Examples
julia> v = Sym[1 0 0; 0 1 0; 0 1 1] ; ℬ = Basis(v)
Basis{3, Sym}
# basis: 3×3 Tensor{2, 3, Sym, 9}:
1 0 0
0 1 0
0 1 1
# dual basis: 3×3 Tensor{2, 3, Sym, 9}:
1 0 0
0 1 -1
0 0 1
# covariant metric tensor: 3×3 SymmetricTensor{2, 3, Sym, 6}:
1 0 0
0 2 1
0 1 1
# contravariant metric tensor: 3×3 SymmetricTensor{2, 3, Sym, 6}:
1 0 0
0 1 -1
0 -1 2
julia> θ, ϕ, ψ = symbols("θ, ϕ, ψ", real = true) ; ℬʳ = Basis(θ, ϕ, ψ) ; display(vecbasis(ℬʳ, :cov))
3×3 Tensor{2, 3, Sym, 9}:
-sin(ψ)⋅sin(ϕ) + cos(θ)⋅cos(ψ)⋅cos(ϕ) -sin(ψ)⋅cos(θ)⋅cos(ϕ) - sin(ϕ)⋅cos(ψ) sin(θ)⋅cos(ϕ)
sin(ψ)⋅cos(ϕ) + sin(ϕ)⋅cos(θ)⋅cos(ψ) -sin(ψ)⋅sin(ϕ)⋅cos(θ) + cos(ψ)⋅cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ)⋅cos(ψ) sin(θ)⋅sin(ψ) cos(θ)
TensND.CanonicalBasis
— TypeCanonicalBasis{dim, T}
Canonical basis of dimension dim
(default: 3) and type T
(default: Sym)
The attributes of this object can be obtained by
vecbasis(ℬ, :cov)
: square matrix defining the primal basiseᵢ=e[:,i]=δᵢⱼ
vecbasis(ℬ, :cont)
: square matrix defining the dual basiseⁱ=E[:,i]=δᵢⱼ
metric(ℬ, :cov)
: square matrix defining the covariant components of the metric tensorgᵢⱼ=eᵢ⋅eⱼ=g[i,j]=δᵢⱼ
metric(ℬ, :cont)
: square matrix defining the contravariant components of the metric tensorgⁱʲ=eⁱ⋅eʲ=gⁱʲ[i,j]=δᵢⱼ
Examples
julia> ℬ = CanonicalBasis()
CanonicalBasis{3, Sym}
# basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# dual basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# covariant metric tensor: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# contravariant metric tensor: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
julia> ℬ₂ = CanonicalBasis{2, Float64}()
CanonicalBasis{2, Float64}
# basis: 2×2 TensND.LazyIdentity{2, Float64}:
1.0 0.0
0.0 1.0
# dual basis: 2×2 TensND.LazyIdentity{2, Float64}:
1.0 0.0
0.0 1.0
# covariant metric tensor: 2×2 TensND.LazyIdentity{2, Float64}:
1.0 0.0
0.0 1.0
# contravariant metric tensor: 2×2 TensND.LazyIdentity{2, Float64}:
1.0 0.0
0.0 1.0
TensND.CoorSystemSym
— TypeCoorSystemSym(OM::AbstractTens{1,dim,Sym},coords::NTuple{dim,Sym},bnorm::AbstractBasis{dim,Sym},χᵢ::NTuple{dim},
tmp_coords::NTuple = (),params::NTuple = ();rules::Dict = Dict(),tmp_var::Dict = Dict(),to_coords::Dict = Dict()) where {dim}
CoorSystemSym(OM::AbstractTens{1,dim,Sym},coords::NTuple{dim,Sym},
tmp_coords::NTuple = (),params::NTuple = ();rules::Dict = Dict(),tmp_var::Dict = Dict(),to_coords::Dict = Dict()) where {dim}
Define a new coordinate system either from
the position vector
OM
, the coordinatescoords
, the basis of unit vectors (𝐞ᵢ
)bnorm
and the Lamé coefficientsχᵢ
In this case the natural basis is formed by the vectors
𝐚ᵢ = χᵢ 𝐞ᵢ
directly calculated from the input data.or the position vector
OM
and the coordinatescoords
In this case the natural basis is formed by the vectors
𝐚ᵢ = ∂ᵢOM
i.e. by the derivative of the position vector with respect to theiᵗʰ
coordinate
Optional parameters can be provided:
tmp_coords
contains temporary variables depending on coordinates (in order to allow symbolic simplifications)params
contains possible parameters involved inOM
rules
contains aDict
with substitution rules to facilitate the simplification of formulastmp_var
contains aDict
with substitution of coordinates by temporary variablesto_coords
indicates how to eliminate the temporary variables to come back to the actual coordinates before derivation for Examples
Examples
julia> ϕ, p = symbols("ϕ p", real = true) ;
julia> p̄, q, q̄, c = symbols("p̄ q q̄ c", positive = true) ;
julia> coords = (ϕ, p, q) ; tmp_coords = (p̄, q̄) ; params = (c,) ;
julia> OM = Tens(c * [p̄ * q̄ * cos(ϕ), p̄ * q̄ * sin(ϕ), p * q]) ;
julia> Spheroidal = CoorSystemSym(OM, coords, tmp_coords, params; tmp_var = Dict(1-p^2 => p̄^2, q^2-1 => q̄^2), to_coords = Dict(p̄ => √(1-p^2), q̄ => √(q^2-1))) ;
TensND.RotatedBasis
— TypeRotatedBasis(θ::T<:Number, ϕ::T<:Number, ψ::T<:Number)
RotatedBasis(θ::T<:Number)
Orthonormal basis of dimension dim
(default: 3) and type T
(default: Sym) built from Euler angles θ
in 2D and (θ, ϕ, ψ)
in 3D
Examples
julia> θ, ϕ, ψ = symbols("θ, ϕ, ψ", real = true) ; ℬʳ = RotatedBasis(θ, ϕ, ψ) ; display(vecbasis(ℬʳ, :cov))
3×3 Tensor{2, 3, Sym, 9}:
-sin(ψ)⋅sin(ϕ) + cos(θ)⋅cos(ψ)⋅cos(ϕ) -sin(ψ)⋅cos(θ)⋅cos(ϕ) - sin(ϕ)⋅cos(ψ) sin(θ)⋅cos(ϕ)
sin(ψ)⋅cos(ϕ) + sin(ϕ)⋅cos(θ)⋅cos(ψ) -sin(ψ)⋅sin(ϕ)⋅cos(θ) + cos(ψ)⋅cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ)⋅cos(ψ) sin(θ)⋅sin(ψ) cos(θ)
TensND.SubManifoldSym
— TypeCoorSystemSym(OM::AbstractTens{1,dim,Sym},coords::NTuple{dim,Sym},bnorm::AbstractBasis{dim,Sym},χᵢ::NTuple{dim},
tmp_coords::NTuple = (),params::NTuple = ();rules::Dict = Dict(),tmp_var::Dict = Dict(),to_coords::Dict = Dict()) where {dim}
CoorSystemSym(OM::AbstractTens{1,dim,Sym},coords::NTuple{dim,Sym},
tmp_coords::NTuple = (),params::NTuple = ();rules::Dict = Dict(),tmp_var::Dict = Dict(),to_coords::Dict = Dict()) where {dim}
Define a new coordinate system either from
the position vector
OM
, the coordinatescoords
, the basis of unit vectors (𝐞ᵢ
)bnorm
and the Lamé coefficientsχᵢ
In this case the natural basis is formed by the vectors
𝐚ᵢ = χᵢ 𝐞ᵢ
directly calculated from the input data.or the position vector
OM
and the coordinatescoords
In this case the natural basis is formed by the vectors
𝐚ᵢ = ∂ᵢOM
i.e. by the derivative of the position vector with respect to theiᵗʰ
coordinate
Optional parameters can be provided:
tmp_coords
contains temporary variables depending on coordinates (in order to allow symbolic simplifications)params
contains possible parameters involved inOM
rules
contains aDict
with substitution rules to facilitate the simplification of formulastmp_var
contains aDict
with substitution of coordinates by temporary variablesto_coords
indicates how to eliminate the temporary variables to come back to the actual coordinates before derivation for Examples
Examples
julia> ϕ, p = symbols("ϕ p", real = true) ;
julia> p̄, q, q̄, c = symbols("p̄ q q̄ c", positive = true) ;
julia> coords = (ϕ, p, q) ; tmp_coords = (p̄, q̄) ; params = (c,) ;
julia> OM = Tens(c * [p̄ * q̄ * cos(ϕ), p̄ * q̄ * sin(ϕ), p * q]) ;
julia> Spheroidal = CoorSystemSym(OM, coords, tmp_coords, params; tmp_var = Dict(1-p^2 => p̄^2, q^2-1 => q̄^2), to_coords = Dict(p̄ => √(1-p^2), q̄ => √(q^2-1))) ;
TensND.Tens
— TypeTens{order,dim,T,A<:AbstractArray}
Tensor type of any order defined by
- a multidata of components (of any type heriting from
AbstractArray
, e.g.Tensor
orSymmetricTensor
) - a basis of
AbstractBasis
type - a tuple of variances (covariant
:cov
or contravariant:cont
) of length equal to theorder
of the tensor
Examples
julia> ℬ = Basis(Sym[1 0 0; 0 1 0; 0 1 1]) ;
julia> T = Tens(metric(ℬ,:cov),ℬ,(:cov,:cov))
Tens{2, 3, Sym, SymmetricTensor{2, 3, Sym, 6}}
# data: 3×3 SymmetricTensor{2, 3, Sym, 6}:
1 0 0
0 2 1
0 1 1
# basis: 3×3 Tensor{2, 3, Sym, 9}:
1 0 0
0 1 0
0 1 1
# var: (:cov, :cov)
julia> components(T,(:cont,:cov),b)
3×3 Matrix{Sym}:
1 0 0
0 1 0
0 0 1
LinearAlgebra.dot
— Methoddot(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a contracted product between two tensors
a ⋅ b = aⁱbⱼ
LinearAlgebra.normalize
— Functionnormalize(ℬ::AbstractBasis, var = cov)
Build a basis after normalization of column vectors of input matrix v
where columns define either
- primal vectors ie
eᵢ=v[:,i]/norm(v[:,i])
ifvar = :cov
as by default - dual vector ie
eⁱ=v[:,i]/norm(v[:,i])
ifvar = :cont
.
TensND.DIV
— MethodDIV(t::AbstractTens{order,dim,Sym},CS::CoorSystemSym{dim}) where {order,dim,T<:Number}
Calculate the divergence of T
with respect to the coordinate system CS
TensND.DIV
— MethodDIV(T::AbstractTens{order,dim,Sym},SM::SubManifoldSym{dim}) where {order,dim}
Calculate the divergence of T
with respect to the coordinate system SM
TensND.GRAD
— MethodGRAD(t::Union{t,AbstractTens{order,dim,T}},CS::CoorSystemSym{dim}) where {order,dim,T<:Number}
Calculate the gradient of t
with respect to the coordinate system CS
TensND.GRAD
— MethodGRAD(T::Union{Sym,AbstractTens{order,dim,Sym}},SM::SubManifoldSym{dim}) where {order,dim}
Calculate the gradient of T
with respect to the coordinate system SM
TensND.HESS
— MethodHESS(t::Union{Sym,AbstractTens{order,dim,Sym}},CS::CoorSystemSym{dim}) where {order,dim,T<:Number}
Calculate the Hessian of T
with respect to the coordinate system CS
TensND.HESS
— MethodHESS(T::Union{Sym,AbstractTens{order,dim,Sym}},SM::SubManifoldSym{dim}) where {order,dim}
Calculate the Hessian of T
with respect to the coordinate system SM
TensND.ISO
— MethodISO(::Val{dim} = Val(3), ::Val{T} = Val(Sym))
Return the three fourth-order isotropic tensors 𝕀, 𝕁, 𝕂
Examples
julia> 𝕀, 𝕁, 𝕂 = ISO() ;
TensND.KM
— MethodKM(t::AbstractTens{order,dim}; kwargs...)
KM(t::AbstractTens{order,dim}, var::NTuple{order,Symbol}, b::AbstractBasis{dim}; kwargs...)
Write the components of a second or fourth order tensor in Kelvin-Mandel notation
Examples
julia> σ = Tens(SymmetricTensor{2,3}((i, j) -> symbols("σ$i$j", real = true))) ;
julia> KM(σ)
6-element Vector{Sym}:
σ11
σ22
σ33
√2⋅σ32
√2⋅σ31
√2⋅σ21
julia> C = Tens(SymmetricTensor{4,3}((i, j, k, l) -> symbols("C$i$j$k$l", real = true))) ;
julia> KM(C)
6×6 Matrix{Sym}:
C₁₁₁₁ C₁₁₂₂ C₁₁₃₃ √2⋅C₁₁₃₂ √2⋅C₁₁₃₁ √2⋅C₁₁₂₁
C₂₂₁₁ C₂₂₂₂ C₂₂₃₃ √2⋅C₂₂₃₂ √2⋅C₂₂₃₁ √2⋅C₂₂₂₁
C₃₃₁₁ C₃₃₂₂ C₃₃₃₃ √2⋅C₃₃₃₂ √2⋅C₃₃₃₁ √2⋅C₃₃₂₁
√2⋅C₃₂₁₁ √2⋅C₃₂₂₂ √2⋅C₃₂₃₃ 2⋅C₃₂₃₂ 2⋅C₃₂₃₁ 2⋅C₃₂₂₁
√2⋅C₃₁₁₁ √2⋅C₃₁₂₂ √2⋅C₃₁₃₃ 2⋅C₃₁₃₂ 2⋅C₃₁₃₁ 2⋅C₃₁₂₁
√2⋅C₂₁₁₁ √2⋅C₂₁₂₂ √2⋅C₂₁₃₃ 2⋅C₂₁₃₂ 2⋅C₂₁₃₁ 2⋅C₂₁₂₁
TensND.KM
— MethodKM(v::AllIsotropic{dim}; kwargs...)
Kelvin-Mandel vector or matrix representation
TensND.LAPLACE
— MethodLAPLACE(t::Union{Sym,AbstractTens{order,dim,Sym}},CS::CoorSystemSym{dim}) where {order,dim,T<:Number}
Calculate the Laplace operator of T
with respect to the coordinate system CS
TensND.LAPLACE
— MethodLAPLACE(T::Union{Sym,AbstractTens{order,dim,Sym}},SM::SubManifoldSym{dim}) where {order,dim}
Calculate the Laplace operator of T
with respect to the coordinate system SM
TensND.LeviCivita
— FunctionLeviCivita(T::Type{<:Number} = Sym)
Builds an Array{T,3} of Levi-Civita Symbol ϵᵢⱼₖ = (i-j) (j-k) (k-i) / 2
Examples
julia> ε = LeviCivita(Sym)
3×3×3 Array{Sym, 3}:
[:, :, 1] =
0 0 0
0 0 1
0 -1 0
[:, :, 2] =
0 0 -1
0 0 0
1 0 0
[:, :, 3] =
0 1 0
-1 0 0
0 0 0
TensND.SYMGRAD
— MethodSYMGRAD(t::Union{T,AbstractTens{order,dim,T}},CS::CoorSystemSym{dim}) where {order,dim,T<:Number}
Calculate the symmetrized gradient of T
with respect to the coordinate system CS
TensND.SYMGRAD
— MethodSYMGRAD(T::Union{Sym,AbstractTens{order,dim,Sym}},SM::SubManifoldSym{dim}) where {order,dim}
Calculate the symmetrized gradient of T
with respect to the coordinate system SM
TensND.angles
— Methodangles(M::AbstractMatrix{T})
Determine the Euler angles corresponding to the input matrix supposed to be a rotation matrix or at least a similarity
Examples
julia> θ, ϕ, ψ = symbols("θ, ϕ, ψ", real = true) ; ℬʳ = RotatedBasis(θ, ϕ, ψ) ; display(vecbasis(ℬʳ, :cov))
3×3 Tensor{2, 3, Sym, 9}:
-sin(ψ)⋅sin(ϕ) + cos(θ)⋅cos(ψ)⋅cos(ϕ) -sin(ψ)⋅cos(θ)⋅cos(ϕ) - sin(ϕ)⋅cos(ψ) sin(θ)⋅cos(ϕ)
sin(ψ)⋅cos(ϕ) + sin(ϕ)⋅cos(θ)⋅cos(ψ) -sin(ψ)⋅sin(ϕ)⋅cos(θ) + cos(ψ)⋅cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ)⋅cos(ψ) sin(θ)⋅sin(ψ) cos(θ)
julia> angles(ℬʳ)
(θ = θ, ϕ = ϕ, ψ = ψ)
TensND.change_tens
— Methodchange_tens(t::AbstractTens{order,dim,T},ℬ::AbstractBasis{dim},var::NTuple{order,Symbol})
change_tens(t::AbstractTens{order,dim,T},ℬ::AbstractBasis{dim})
change_tens(t::AbstractTens{order,dim,T},var::NTuple{order,Symbol})
Rewrite the same tensor with components corresponding to new variances and/or to a new basis
julia> ℬ = Basis(Sym[0 1 1; 1 0 1; 1 1 0]) ;
julia> TV = Tens(Tensor{1,3}(i->symbols("v$i",real=true)))
TensND.TensCanonical{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
v₁
v₂
v₃
# basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont,)
julia> factor.(components(TV, ℬ, (:cont,)))
3-element Vector{Sym}:
-(v1 - v2 - v3)/2
(v1 - v2 + v3)/2
(v1 + v2 - v3)/2
julia> ℬ₀ = Basis(Sym[0 1 1; 1 0 1; 1 1 1]) ;
julia> TV0 = change_tens(TV, ℬ₀)
Tens{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
-v₁ + v₃
-v₂ + v₃
v₁ + v₂ - v₃
# basis: 3×3 Tensor{2, 3, Sym, 9}:
0 1 1
1 0 1
1 1 1
# var: (:cont,)
TensND.change_tens_canon
— Methodchange_tens_canon(t::AbstractTens{order,dim,T},var::NTuple{order,Symbol})
Rewrite the same tensor with components corresponding to the canonical basis
julia> ℬ = Basis(Sym[0 1 1; 1 0 1; 1 1 0]) ;
julia> TV = Tens(Tensor{1,3}(i->symbols("v$i",real=true)), ℬ)
Tens{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
v₁
v₂
v₃
# basis: 3×3 Tensor{2, 3, Sym, 9}:
0 1 1
1 0 1
1 1 1
# var: (:cont,)
julia> TV0 = change_tens_canon(TV)
TensND.TensCanonical{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
v₂ + v₃
v₁ + v₃
v₁ + v₂ + v₃
# basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont,)
TensND.components
— Methodcomponents(t::AbstractTens{order,dim,T},ℬ::AbstractBasis{dim},var::NTuple{order,Symbol})
components(t::AbstractTens{order,dim,T},ℬ::AbstractBasis{dim})
components(t::AbstractTens{order,dim,T},var::NTuple{order,Symbol})
Extract the components of a tensor for new variances and/or in a new basis
Examples
julia> ℬ = Basis(Sym[0 1 1; 1 0 1; 1 1 0]) ;
julia> TV = Tens(Tensor{1,3}(i->symbols("v$i",real=true)))
TensND.TensCanonical{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
v₁
v₂
v₃
# basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont,)
julia> factor.(components(TV, ℬ, (:cont,)))
3-element Vector{Sym}:
-(v1 - v2 - v3)/2
(v1 - v2 + v3)/2
(v1 + v2 - v3)/2
julia> components(TV, ℬ, (:cov,))
3-element Vector{Sym}:
v₂ + v₃
v₁ + v₃
v₁ + v₂
julia> simplify.(components(TV, normalize(ℬ), (:cov,)))
3-element Vector{Sym}:
sqrt(2)*(v2 + v3)/2
sqrt(2)*(v1 + v3)/2
sqrt(2)*(v1 + v2)/2
julia> TT = Tens(Tensor{2,3}((i,j)->symbols("t$i$j",real=true)))
TensND.TensCanonical{2, 3, Sym, Tensor{2, 3, Sym, 9}}
# data: 3×3 Tensor{2, 3, Sym, 9}:
t₁₁ t₁₂ t₁₃
t₂₁ t₂₂ t₂₃
t₃₁ t₃₂ t₃₃
# basis: 3×3 TensND.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont, :cont)
julia> components(TT, ℬ, (:cov,:cov))
3×3 Matrix{Sym}:
t₂₂ + t₂₃ + t₃₂ + t₃₃ t₂₁ + t₂₃ + t₃₁ + t₃₃ t₂₁ + t₂₂ + t₃₁ + t₃₂
t₁₂ + t₁₃ + t₃₂ + t₃₃ t₁₁ + t₁₃ + t₃₁ + t₃₃ t₁₁ + t₁₂ + t₃₁ + t₃₂
t₁₂ + t₁₃ + t₂₂ + t₂₃ t₁₁ + t₁₃ + t₂₁ + t₂₃ t₁₁ + t₁₂ + t₂₁ + t₂₂
julia> factor.(components(TT, ℬ, (:cont,:cov)))
3×3 Matrix{Sym}:
-(t12 + t13 - t22 - t23 - t32 - t33)/2 … -(t11 + t12 - t21 - t22 - t31 - t32)/2
(t12 + t13 - t22 - t23 + t32 + t33)/2 (t11 + t12 - t21 - t22 + t31 + t32)/2
(t12 + t13 + t22 + t23 - t32 - t33)/2 (t11 + t12 + t21 + t22 - t31 - t32)/2
TensND.components_canon
— Methodcomponents_canon(t::AbstractTens)
Extract the components of a tensor in the canonical basis
TensND.contract
— Methodcontract(t::AbstractTens{order,dim}, i::Integer, j::Integer)
Calculate the tensor obtained after contraction with respect to the indices i
and j
TensND.coorsys_cartesian
— Methodcoorsys_cartesian(coords = symbols("x y z", real = true))
Return the cartesian coordinate system
Examples
julia> Cartesian = coorsys_cartesian() ; 𝐗 = getcoords(Cartesian) ; 𝐄 = unitvec(Cartesian) ; ℬ = getbasis(Cartesian)
julia> 𝛔 = Tens(SymmetricTensor{2,3}((i, j) -> SymFunction("σ$i$j", real = true)(𝐗...))) ;
julia> DIV(𝛔, CScar)
Tens.TensCanonical{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
Derivative(σ11(x, y, z), x) + Derivative(σ21(x, y, z), y) + Derivative(σ31(x, y, z), z)
Derivative(σ21(x, y, z), x) + Derivative(σ22(x, y, z), y) + Derivative(σ32(x, y, z), z)
Derivative(σ31(x, y, z), x) + Derivative(σ32(x, y, z), y) + Derivative(σ33(x, y, z), z)
# basis: 3×3 Tens.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont,)
TensND.coorsys_cylindrical
— Methodcoorsys_cylindrical(coords = (symbols("r", positive = true), symbols("θ", real = true), symbols("z", real = true)); canonical = false)
Return the cylindrical coordinate system
Examples
julia> Cylindrical = coorsys_cylindrical() ; rθz = getcoords(Cylindrical) ; 𝐞ʳ, 𝐞ᶿ, 𝐞ᶻ = unitvec(Cylindrical) ; ℬᶜ = getbasis(Cylindrical)
julia> 𝐯 = Tens(Vec{3}(i -> SymFunction("v$(rθz[i])", real = true)(rθz...)), ℬᶜ) ;
julia> DIV(𝐯, Cylindrical)
∂
vr(r, θ, z) + ──(vθ(r, θ, z))
∂ ∂ ∂θ
──(vr(r, θ, z)) + ──(vz(r, θ, z)) + ─────────────────────────────
∂r ∂z r
TensND.coorsys_polar
— Methodcoorsys_polar(coords = (symbols("r", positive = true), symbols("θ", real = true)); canonical = false)
Return the polar coordinate system
Examples
julia> Polar = coorsys_polar() ; r, θ = getcoords(Polar) ; 𝐞ʳ, 𝐞ᶿ = unitvec(Polar) ; ℬᵖ = getbasis(Polar)
julia> f = SymFunction("f", real = true)(r, θ) ;
julia> LAPLACE(f, Polar)
2
∂
───(f(r, θ))
2
∂ ∂θ
2 ──(f(r, θ)) + ────────────
∂ ∂r r
───(f(r, θ)) + ──────────────────────────
2 r
∂r
TensND.coorsys_spherical
— Methodcoorsys_spherical(coords = (symbols("θ", real = true), symbols("ϕ", real = true), symbols("r", positive = true)); canonical = false)
Return the spherical coordinate system
Examples
julia> Spherical = coorsys_spherical() ; θ, ϕ, r = getcoords(Spherical) ; 𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ = unitvec(Spherical) ; ℬˢ = getbasis(Spherical)
julia> for σⁱʲ ∈ ("σʳʳ", "σᶿᶿ", "σᵠᵠ") @eval $(Symbol(σⁱʲ)) = SymFunction($σⁱʲ, real = true)($r) end ;
julia> 𝛔 = σʳʳ * 𝐞ʳ ⊗ 𝐞ʳ + σᶿᶿ * 𝐞ᶿ ⊗ 𝐞ᶿ + σᵠᵠ * 𝐞ᵠ ⊗ 𝐞ᵠ ;
julia> div𝛔 = DIV(𝛔, Spherical)
Tens.TensRotated{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
(-σᵠᵠ(r) + σᶿᶿ(r))*cos(θ)/(r*sin(θ))
0
Derivative(σʳʳ(r), r) + (σʳʳ(r) - σᵠᵠ(r))/r + (σʳʳ(r) - σᶿᶿ(r))/r
# basis: 3×3 Tensor{2, 3, Sym, 9}:
cos(θ)⋅cos(ϕ) -sin(ϕ) sin(θ)⋅cos(ϕ)
sin(ϕ)⋅cos(θ) cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ) 0 cos(θ)
# var: (:cont,)
julia> div𝛔 ⋅ 𝐞ʳ
d σʳʳ(r) - σᵠᵠ(r) σʳʳ(r) - σᶿᶿ(r)
──(σʳʳ(r)) + ─────────────── + ───────────────
dr r r
TensND.coorsys_spheroidal
— Methodcoorsys_spheroidal(coords = (symbols("ϕ", real = true),symbols("p", real = true),symbols("q", positive = true),),
c = symbols("c", positive = true),tmp_coords = (symbols("p̄ q̄", positive = true)...,),)
Return the spheroidal coordinate system
Examples
julia> Spheroidal = coorsys_spheroidal() ; OM = getOM(Spheroidal)
Tens.TensCanonical{1, 3, Sym, Vec{3, Sym}}
# data: 3-element Vec{3, Sym}:
c⋅p̄⋅q̄⋅cos(ϕ)
c⋅p̄⋅q̄⋅sin(ϕ)
c⋅p⋅q
# basis: 3×3 Tens.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
# var: (:cont,)
julia> LAPLACE(OM[1]^2, Spheroidal)
2
TensND.init_cartesian
— Functioninit_cartesian(coords = symbols("x y z", real = true))
Returns the coordinates, unit vectors and basis of the cartesian basis
Examples
julia> coords, vectors, ℬ = init_cartesian() ; x, y, z = coords ; 𝐞₁, 𝐞₂, 𝐞₃ = vectors ;
TensND.init_cylindrical
— Functioninit_cylindrical(coords = (symbols("r", positive = true), symbols("θ z", real = true)...); canonical = false)
Returns the coordinates, base vectors and basis of the cylindrical basis
Examples
julia> coords, vectors, ℬᶜ = init_cylindrical() ; r, θ, z = coords ; 𝐞ʳ, 𝐞ᶿ, 𝐞ᶻ = vectors ;
TensND.init_polar
— Functioninit_polar(coords = (symbols("r θ", real = true)); canonical = false)
Returns the coordinates, base vectors and basis of the polar basis
Examples
julia> coords, vectors, ℬᵖ = init_polar() ; r, θ = coords ; 𝐞ʳ, 𝐞ᶿ = vectors ;
TensND.init_rotated
— Functioninit_rotated(coords = symbols("θ ϕ ψ", real = true); canonical = false)
Return the angles, base vectors and basis of the rotated basis. Note that here the coordinates are angles and do not represent a valid parametrization of ℝ³
Examples
julia> angles, vectors, ℬʳ = init_rotated() ; θ, ϕ, ψ = angles ; 𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ = vectors ;
TensND.init_spherical
— Functioninit_spherical(coords = (symbols("θ ϕ", real = true)..., symbols("r", positive = true)); canonical = false)
Return the coordinates, base vectors and basis of the spherical basis. Take care that the order of the 3 vectors is 𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ
so that the basis coincides with the canonical one when the angles are null and in consistency the coordinates are ordered as θ, ϕ, r
.
Examples
julia> coords, vectors, ℬˢ = init_spherical() ; θ, ϕ, r = coords ; 𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ = vectors ;
TensND.invKM
— MethodinvKM(v::AbstractVecOrMat; kwargs...)
Define a tensor from a Kelvin-Mandel vector or matrix representation
TensND.isorthogonal
— Methodisorthogonal(ℬ::AbstractBasis)
Check whether the basis ℬ
is orthogonal
TensND.isorthonormal
— Methodisorthonormal(ℬ::AbstractBasis)
Check whether the basis ℬ
is orthonormal
TensND.metric
— Methodmetric(ℬ::AbstractBasis, var = :cov)
Return the covariant (if var = :cov
) or contravariant (if var = :cont
) metric matrix
TensND.otimesul
— Methodotimesul(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a special tensor product between two tensors of at least second order
(𝐚 ⊠ˢ 𝐛) ⊡ 𝐩 = (𝐚 ⊠ 𝐛) ⊡ (𝐩 + ᵗ𝐩)/2 = 1/2(aⁱᵏbʲˡ+aⁱˡbʲᵏ) pₖₗ eᵢ⊗eⱼ
TensND.qcontract
— Methodqcontract(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a quadruple contracted product between two tensors
𝔸 ⊙ 𝔹 = AᵢⱼₖₗBⁱʲᵏˡ
Examples
julia> 𝕀 = t𝕀(Sym) ; 𝕁 = t𝕁(Sym) ; 𝕂 = t𝕂(Sym) ;
julia> 𝕀 ⊙ 𝕀
6
julia> 𝕁 ⊙ 𝕀
1
julia> 𝕂 ⊙ 𝕀
5
julia> 𝕂 ⊙ 𝕁
0
TensND.rot2
— Methodrot2(θ)
Return a 2D rotation matrix with respect to the angle θ
Examples
julia> rot2(θ)
2×2 Tensor{2, 2, Sym, 4}:
cos(θ) -sin(θ)
sin(θ) cos(θ)
TensND.rot3
— Functionrot3(θ, ϕ = 0, ψ = 0)
Return a rotation matrix with respect to the 3 Euler angles θ, ϕ, ψ
Examples
julia> cθ, cϕ, cψ, sθ, sϕ, sψ = symbols("cθ cϕ cψ sθ sϕ sψ", real = true) ;
julia> d = Dict(cos(θ) => cθ, cos(ϕ) => cϕ, cos(ψ) => cψ, sin(θ) => sθ, sin(ϕ) => sϕ, sin(ψ) => sψ) ;
julia> subs.(rot3(θ, ϕ, ψ),d...)
3×3 StaticArrays.SMatrix{3, 3, Sym, 9} with indices SOneTo(3)×SOneTo(3):
cθ⋅cψ⋅cϕ - sψ⋅sϕ -cθ⋅cϕ⋅sψ - cψ⋅sϕ cϕ⋅sθ
cθ⋅cψ⋅sϕ + cϕ⋅sψ -cθ⋅sψ⋅sϕ + cψ⋅cϕ sθ⋅sϕ
-cψ⋅sθ sθ⋅sψ cθ
TensND.rot6
— Functionrot6(θ, ϕ = 0, ψ = 0)
Return a rotation matrix with respect to the 3 Euler angles θ, ϕ, ψ
Examples
julia> cθ, cϕ, cψ, sθ, sϕ, sψ = symbols("cθ cϕ cψ sθ sϕ sψ", real = true) ;
julia> d = Dict(cos(θ) => cθ, cos(ϕ) => cϕ, cos(ψ) => cψ, sin(θ) => sθ, sin(ϕ) => sϕ, sin(ψ) => sψ) ;
julia> R = Tens(subs.(rot3(θ, ϕ, ψ),d...))
Tens.TensCanonical{2, 3, Sym, Tensor{2, 3, Sym, 9}}
# data: 3×3 Tensor{2, 3, Sym, 9}:
cθ⋅cψ⋅cϕ - sψ⋅sϕ -cθ⋅cϕ⋅sψ - cψ⋅sϕ cϕ⋅sθ
cθ⋅cψ⋅sϕ + cϕ⋅sψ -cθ⋅sψ⋅sϕ + cψ⋅cϕ sθ⋅sϕ
-cψ⋅sθ sθ⋅sψ cθ
# var: (:cont, :cont)
# basis: 3×3 Tens.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
julia> RR = R ⊠ˢ R
Tens.TensCanonical{4, 3, Sym, SymmetricTensor{4, 3, Sym, 36}}
# data: 6×6 Matrix{Sym}:
(cθ*cψ*cϕ - sψ*sϕ)^2 (-cθ*cϕ*sψ - cψ*sϕ)^2 cϕ^2*sθ^2 √2⋅cϕ⋅sθ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) √2⋅cϕ⋅sθ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) √2⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ)
(cθ*cψ*sϕ + cϕ*sψ)^2 (-cθ*sψ*sϕ + cψ*cϕ)^2 sθ^2*sϕ^2 √2⋅sθ⋅sϕ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) √2⋅sθ⋅sϕ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ)⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ)
cψ^2*sθ^2 sθ^2*sψ^2 cθ^2 √2⋅cθ⋅sθ⋅sψ -√2⋅cθ⋅cψ⋅sθ -sqrt(2)*cψ*sθ^2*sψ
-√2⋅cψ⋅sθ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅sθ⋅sψ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) √2⋅cθ⋅sθ⋅sϕ cθ*(-cθ*sψ*sϕ + cψ*cϕ) + sθ^2*sψ*sϕ cθ*(cθ*cψ*sϕ + cϕ*sψ) - cψ*sθ^2*sϕ -cψ⋅sθ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) + sθ⋅sψ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ)
-√2⋅cψ⋅sθ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) √2⋅sθ⋅sψ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) √2⋅cθ⋅cϕ⋅sθ cθ*(-cθ*cϕ*sψ - cψ*sϕ) + cϕ*sθ^2*sψ cθ*(cθ*cψ*cϕ - sψ*sϕ) - cψ*cϕ*sθ^2 -cψ⋅sθ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) + sθ⋅sψ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)
√2⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ)⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) sqrt(2)*cϕ*sθ^2*sϕ cϕ⋅sθ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) + sθ⋅sϕ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) cϕ⋅sθ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) + sθ⋅sϕ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) (cθ*cψ*cϕ - sψ*sϕ)*(-cθ*sψ*sϕ + cψ*cϕ) + (cθ*cψ*sϕ + cϕ*sψ)*(-cθ*cϕ*sψ - cψ*sϕ)
# var: (:cont, :cont, :cont, :cont)
# basis: 3×3 Tens.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
julia> R6 = invKM(subs.(KM(rot6(θ, ϕ, ψ)),d...))
Tens.TensCanonical{4, 3, Sym, SymmetricTensor{4, 3, Sym, 36}}
# data: 6×6 Matrix{Sym}:
(cθ*cψ*cϕ - sψ*sϕ)^2 (-cθ*cϕ*sψ - cψ*sϕ)^2 cϕ^2*sθ^2 √2⋅cϕ⋅sθ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) √2⋅cϕ⋅sθ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) √2⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ)
(cθ*cψ*sϕ + cϕ*sψ)^2 (-cθ*sψ*sϕ + cψ*cϕ)^2 sθ^2*sϕ^2 √2⋅sθ⋅sϕ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) √2⋅sθ⋅sϕ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ)⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ)
cψ^2*sθ^2 sθ^2*sψ^2 cθ^2 √2⋅cθ⋅sθ⋅sψ -√2⋅cθ⋅cψ⋅sθ -sqrt(2)*cψ*sθ^2*sψ
-√2⋅cψ⋅sθ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅sθ⋅sψ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) √2⋅cθ⋅sθ⋅sϕ cθ*(-cθ*sψ*sϕ + cψ*cϕ) + sθ^2*sψ*sϕ cθ*(cθ*cψ*sϕ + cϕ*sψ) - cψ*sθ^2*sϕ -cψ⋅sθ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) + sθ⋅sψ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ)
-√2⋅cψ⋅sθ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) √2⋅sθ⋅sψ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) √2⋅cθ⋅cϕ⋅sθ cθ*(-cθ*cϕ*sψ - cψ*sϕ) + cϕ*sθ^2*sψ cθ*(cθ*cψ*cϕ - sψ*sϕ) - cψ*cϕ*sθ^2 -cψ⋅sθ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) + sθ⋅sψ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)
√2⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ)⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) √2⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ)⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) sqrt(2)*cde Liv Lehn ϕ*sθ^2*sϕ cϕ⋅sθ⋅(-cθ⋅sψ⋅sϕ + cψ⋅cϕ) + sθ⋅sϕ⋅(-cθ⋅cϕ⋅sψ - cψ⋅sϕ) cϕ⋅sθ⋅(cθ⋅cψ⋅sϕ + cϕ⋅sψ) + sθ⋅sϕ⋅(cθ⋅cψ⋅cϕ - sψ⋅sϕ) (cθ*cψ*cϕ - sψ*sϕ)*(-cθ*sψ*sϕ + cψ*cϕ) + (cθ*cψ*sϕ + cϕ*sψ)*(-cθ*cϕ*sψ - cψ*sϕ)
# var: (:cont, :cont, :cont, :cont)
# basis: 3×3 Tens.LazyIdentity{3, Sym}:
1 0 0
0 1 0
0 0 1
julia> R6 == RR
true
TensND.sotimes
— Methodsotimes(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a symmetric tensor product between two tensors
(aⁱeᵢ) ⊗ˢ (bʲeⱼ) = 1/2(aⁱbʲ + aʲbⁱ) eᵢ⊗eⱼ
TensND.tensId2
— MethodtensId2(::Val{dim}, ::Val{T}) where {dim,T<:Number}
Identity tensor of second order 𝟏ᵢⱼ = δᵢⱼ = 1 if i=j otherwise 0
Examples
julia> 𝟏 = t𝟏() ; KM(𝟏)
6-element Vector{Sym}:
1
1
1
0
0
0
julia> 𝟏.data
3×3 SymmetricTensor{2, 3, Sym, 6}:
1 0 0
0 1 0
0 0 1
TensND.tensId4
— MethodtensId4(::Val{dim} = Val(3), ::Val{T} = Val(Sym))
Symmetric identity tensor of fourth order 𝕀 = 𝟏 ⊠ˢ 𝟏
i.e. (𝕀)ᵢⱼₖₗ = (δᵢₖδⱼₗ+δᵢₗδⱼₖ)/2
Examples
julia> 𝕀 = t𝕀() ; KM(𝕀)
6×6 Matrix{Sym}:
1 0 0 0 0 0
0 1 0 0 0 0
0 0 1 0 0 0
0 0 0 1 0 0
0 0 0 0 1 0
0 0 0 0 0 1
TensND.tensJ4
— MethodtensJ4(::Val{dim} = Val(3), ::Val{T} = Val(Sym))
Spherical projector of fourth order 𝕁 = (𝟏 ⊗ 𝟏) / dim
i.e. (𝕁)ᵢⱼₖₗ = δᵢⱼδₖₗ/dim
Examples
julia> 𝕁 = t𝕁() ; KM(𝕁)
6×6 Matrix{Sym}:
1/3 1/3 1/3 0 0 0
1/3 1/3 1/3 0 0 0
1/3 1/3 1/3 0 0 0
0 0 0 0 0 0
0 0 0 0 0 0
0 0 0 0 0 0
TensND.tensK4
— MethodtensK4(::Val{dim} = Val(3), ::Val{T} = Val(Sym))
Deviatoric projector of fourth order 𝕂 = 𝕀 - 𝕁
i.e. (𝕂)ᵢⱼₖₗ = (δᵢₖδⱼₗ+δᵢₗδⱼₖ)/2 - δᵢⱼδₖₗ/dim
Examples
julia> 𝕂 = t𝕂() ; KM(𝕂)
6×6 Matrix{Sym}:
2/3 -1/3 -1/3 0 0 0
-1/3 2/3 -1/3 0 0 0
-1/3 -1/3 2/3 0 0 0
0 0 0 1 0 0
0 0 0 0 1 0
0 0 0 0 0 1
TensND.vecbasis
— Methodvecbasis(ℬ::AbstractBasis, var = :cov)
Return the primal (if var = :cov
) or dual (if var = :cont
) basis
TensND.∂
— Method∂(t::AbstractTens{order,dim,T,A},xᵢ::T)
Return the derivative of the tensor t
with respect to the variable x_i
Examples
julia> (θ, ϕ, r), (𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ), ℬˢ = init_spherical() ;
julia> ∂(𝐞ʳ, ϕ) == sin(θ) * 𝐞ᵠ
true
julia> ∂(𝐞ʳ ⊗ 𝐞ʳ,θ)
Tens.TensRotated{2, 3, Sym, SymmetricTensor{2, 3, Sym, 6}}
# data: 3×3 SymmetricTensor{2, 3, Sym, 6}:
0 0 1
0 0 0
1 0 0
# basis: 3×3 Tensor{2, 3, Sym, 9}:
cos(θ)⋅cos(ϕ) -sin(ϕ) sin(θ)⋅cos(ϕ)
sin(ϕ)⋅cos(θ) cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ) 0 cos(θ)
# var: (:cont, :cont)
TensND.𝐞
— Method𝐞(i::Integer, dim::Int = 3, T::Type{<:Number} = Sym)
Vector of the canonical basis
Examples
julia> 𝐞(1)
Tens{1, 3, Sym, Sym, Vec{3, Sym}, CanonicalBasis{3, Sym}}
# data: 3-element Vec{3, Sym}:
1
0
0
# var: (:cont,)
# basis: 3×3 Tensor{2, 3, Sym, 9}:
1 0 0
0 1 0
0 0 1
TensND.𝐞ˢ
— Method𝐞ˢ(i::Integer, θ::T = zero(Sym), ϕ::T = zero(Sym), ψ::T = zero(Sym); canonical = false)
Vector of the basis rotated with the 3 Euler angles θ, ϕ, ψ
(spherical if ψ=0
)
Examples
julia> θ, ϕ, ψ = symbols("θ, ϕ, ψ", real = true) ;
Tens{1, 3, Sym, Sym, Vec{3, Sym}, RotatedBasis{3, Sym}}
# data: 3-element Vec{3, Sym}:
1
0
0
# var: (:cont,)
# basis: 3×3 Tensor{2, 3, Sym, 9}:
-sin(ψ)⋅sin(ϕ) + cos(θ)⋅cos(ψ)⋅cos(ϕ) -sin(ψ)⋅cos(θ)⋅cos(ϕ) - sin(ϕ)⋅cos(ψ) sin(θ)⋅cos(ϕ)
sin(ψ)⋅cos(ϕ) + sin(ϕ)⋅cos(θ)⋅cos(ψ) -sin(ψ)⋅sin(ϕ)⋅cos(θ) + cos(ψ)⋅cos(ϕ) sin(θ)⋅sin(ϕ)
-sin(θ)⋅cos(ψ) sin(θ)⋅sin(ψ) cos(θ)
TensND.𝐞ᵖ
— Method𝐞ᵖ(i::Integer, θ::T = zero(Sym); canonical = false)
Vector of the polar basis
Examples
julia> θ = symbols("θ", real = true) ;
julia> 𝐞ᵖ(1, θ)
Tens{1, 2, Sym, Sym, Vec{2, Sym}, RotatedBasis{2, Sym}}
# data: 2-element Vec{2, Sym}:
1
0
# var: (:cont,)
# basis: 2×2 Tensor{2, 2, Sym, 4}:
cos(θ) -sin(θ)
sin(θ) cos(θ)
TensND.𝐞ᶜ
— Method𝐞ᶜ(i::Integer, θ::T = zero(Sym); canonical = false)
Vector of the cylindrical basis
Examples
julia> θ = symbols("θ", real = true) ;
julia> 𝐞ᶜ(1, θ)
Tens{1, 3, Sym, Sym, Vec{3, Sym}, RotatedBasis{3, Sym}}
# data: 3-element Vec{3, Sym}:
1
0
0
# var: (:cont,)
# basis: 3×3 Tensor{2, 3, Sym, 9}:
cos(θ) -sin(θ) 0
sin(θ) cos(θ) 0
0 0 1
Tensors.dcontract
— Methoddcontract(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a double contracted product between two tensors
𝛔 ⊡ 𝛆 = σⁱʲεᵢⱼ
𝛔 = ℂ ⊡ 𝛆
Examples
julia> 𝛆 = Tens(SymmetricTensor{2,3}((i, j) -> symbols("ε$i$j", real = true))) ;
julia> k, μ = symbols("k μ", real =true) ;
julia> ℂ = 3k * t𝕁() + 2μ * t𝕂() ;
julia> 𝛔 = ℂ ⊡ 𝛆
Tens{2, 3, Sym, Sym, SymmetricTensor{2, 3, Sym, 6}, CanonicalBasis{3, Sym}}
# data: 3×3 SymmetricTensor{2, 3, Sym, 6}:
ε11*(k + 4*μ/3) + ε22*(k - 2*μ/3) + ε33*(k - 2*μ/3) 2⋅ε21⋅μ 2⋅ε31⋅μ
2⋅ε21⋅μ ε11*(k - 2*μ/3) + ε22*(k + 4*μ/3) + ε33*(k - 2*μ/3) 2⋅ε32⋅μ
2⋅ε31⋅μ 2⋅ε32⋅μ ε11*(k - 2*μ/3) + ε22*(k - 2*μ/3) + ε33*(k + 4*μ/3)
# var: (:cont, :cont)
# basis: 3×3 Tensor{2, 3, Sym, 9}:
1 0 0
0 1 0
0 0 1
Tensors.dotdot
— Methoddotdot(v1::AbstractTens{order1,dim}, S::AbstractTens{orderS,dim}, v2::AbstractTens{order2,dim})
Define a bilinear operator 𝐯₁⋅𝕊⋅𝐯₂
Examples
julia> n = Tens(Sym[0, 0, 1]) ;
julia> k, μ = symbols("k μ", real =true) ;
julia> ℂ = 3k * t𝕁() + 2μ * t𝕂() ;
julia> dotdot(n,ℂ,n) # Acoustic tensor
3×3 Tens{2, 3, Sym, Sym, Tensor{2, 3, Sym, 9}, CanonicalBasis{3, Sym}}:
μ 0 0
0 μ 0
0 0 k + 4*μ/3
Tensors.otimes
— Methodotimes(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a tensor product between two tensors
(aⁱeᵢ) ⊗ (bʲeⱼ) = aⁱbʲ eᵢ⊗eⱼ
Tensors.otimesu
— Methodotimesu(t1::AbstractTens{order1,dim}, t2::AbstractTens{order2,dim})
Define a special tensor product between two tensors of at least second order
(𝐚 ⊠ 𝐛) ⊡ 𝐩 = 𝐚⋅𝐩⋅𝐛 = aⁱᵏbʲˡpₖₗ eᵢ⊗eⱼ
TensND.@set_coorsys
— Macro@set_coorsys CS
@set_coorsys(CS)
Set a coordinate system in order to avoid precising it in differential operators
Examples
julia> Spherical = coorsys_spherical() ; θ, ϕ, r = getcoords(Spherical) ; 𝐞ᶿ, 𝐞ᵠ, 𝐞ʳ = unitvec(Spherical) ; vec = ("𝐞ᶿ", "𝐞ᵠ", "𝐞ʳ") ;
julia> @set_coorsys Spherical
julia> intrinsic(GRAD(𝐞ʳ),vec)
(1/r)𝐞ᶿ⊗𝐞ᶿ + (1/r)𝐞ᵠ⊗𝐞ᵠ
julia> intrinsic(DIV(𝐞ʳ ⊗ 𝐞ʳ),vec)
(2/r)𝐞ʳ
julia> LAPLACE(1/r)
0
Index
TensND.Basis
TensND.CanonicalBasis
TensND.CoorSystemSym
TensND.RotatedBasis
TensND.SubManifoldSym
TensND.Tens
LinearAlgebra.dot
LinearAlgebra.normalize
TensND.DIV
TensND.DIV
TensND.GRAD
TensND.GRAD
TensND.HESS
TensND.HESS
TensND.ISO
TensND.KM
TensND.KM
TensND.LAPLACE
TensND.LAPLACE
TensND.LeviCivita
TensND.SYMGRAD
TensND.SYMGRAD
TensND.angles
TensND.change_tens
TensND.change_tens_canon
TensND.components
TensND.components_canon
TensND.contract
TensND.coorsys_cartesian
TensND.coorsys_cylindrical
TensND.coorsys_polar
TensND.coorsys_spherical
TensND.coorsys_spheroidal
TensND.init_cartesian
TensND.init_cylindrical
TensND.init_polar
TensND.init_rotated
TensND.init_spherical
TensND.invKM
TensND.isorthogonal
TensND.isorthonormal
TensND.metric
TensND.otimesul
TensND.qcontract
TensND.rot2
TensND.rot3
TensND.rot6
TensND.sotimes
TensND.tensId2
TensND.tensId4
TensND.tensJ4
TensND.tensK4
TensND.vecbasis
TensND.∂
TensND.𝐞
TensND.𝐞ˢ
TensND.𝐞ᵖ
TensND.𝐞ᶜ
Tensors.dcontract
Tensors.dotdot
Tensors.otimes
Tensors.otimesu
TensND.@set_coorsys