PlasmaEquilibriumToolkit
PlasmaEquilibriumToolkit.AbstractMagneticCoordinates
PlasmaEquilibriumToolkit.AbstractMagneticEquilibrium
PlasmaEquilibriumToolkit.BasisVectors
PlasmaEquilibriumToolkit.BoozerCoordinates
PlasmaEquilibriumToolkit.ClebschCoordinates
PlasmaEquilibriumToolkit.CoordinateVector
PlasmaEquilibriumToolkit.FluxCoordinates
PlasmaEquilibriumToolkit.NullEquilibrium
PlasmaEquilibriumToolkit.PestCoordinates
Base.abs
PlasmaEquilibriumToolkit.MagneticCoordinateGrid
PlasmaEquilibriumToolkit.basis_vectors
PlasmaEquilibriumToolkit.curvature_components
PlasmaEquilibriumToolkit.geodesic_curvature
PlasmaEquilibriumToolkit.grad_B_projection
PlasmaEquilibriumToolkit.jacobian
PlasmaEquilibriumToolkit.metric
PlasmaEquilibriumToolkit.normal_curvature
PlasmaEquilibriumToolkit.transform_basis
PlasmaEquilibriumToolkit.transform_basis
PlasmaEquilibriumToolkit.AbstractMagneticCoordinates
— TypeAbstractMagneticCoordinates
Abstract supertype for different magnetic coordinates.
PlasmaEquilibriumToolkit.AbstractMagneticEquilibrium
— TypeAbstractMagneticEquilibrium
Abstract supertype for different magnetic equilibrium representations (VMEC, SPEC,...).
PlasmaEquilibriumToolkit.BasisVectors
— TypeBasisVectors{T} <: SArray{Tuple{3,3},T,2,9}
Convenience type for defining the Cartesian coordinate representation of 3D curvilinear vector fields. For a vector defined by inverse map $R(x,y,z) = (U(x,y,z),V(x,y,z),W(x,y,z))$, the components are represented by the 3×3 SArray:
R = [Ux Vx Wx;Uy Vy Wy;Uz Vz Wz]
where Ux
is the x
-component of U
.
Examples
julia> using StaticArrays
julia> typeof(@SArray(ones(3,3))) <: BasisVectors{Float64}
true
PlasmaEquilibriumToolkit.BoozerCoordinates
— TypeBoozerCoordinates{T,A}(ψ::T,χ::A,ϕ::A) <: MagneticCoordinates
Coordinates on a magnetic flux surface, where ψ
is the flux label and χ
and ϕ
are angle-like variables
PlasmaEquilibriumToolkit.ClebschCoordinates
— TypeClebschCoordinates(α,β,η)
Coordintes (α,β,η)
representing a divergence-free magnetic field $B = ∇α×∇β$ with a Jacobian defined by $√g = 1/(∇η ⋅ ∇α × ∇β)$.
PlasmaEquilibriumToolkit.CoordinateVector
— TypeCoordinateVector{T} <: SVector{3,T}
Convenience type for defining a 3-element StaticVector.
PlasmaEquilibriumToolkit.FluxCoordinates
— TypeFluxCoordinates{T,A}(ψ::T,θ::A,ζ::A) <: MagneticCoordinates
Coordinates on a magnetic flux surface, where ψ
is the physical toroidal flux divided by 2π and θ
and ζ
are angle-like variables
PlasmaEquilibriumToolkit.NullEquilibrium
— TypeNullEquilibrium()
Empty subtype of MagneticEquilibrium to represent no equilibrium
PlasmaEquilibriumToolkit.PestCoordinates
— TypePestCoordinates{T,A}(ψ::T,α::A,ζ::A) <: MagneticCoordinates
Coordinates on a magnetic flux surface, where ψ
is the toroidal flux divided by 2π in the clockwise direction, α
is the field line label, and ζ
is the geometric toroidal angle advancing the clockwise direction.
Base.abs
— Methodabs(e::BasisVectors{T},component=0)
Compute the L2 norm of the basis vector e
. If component
> 0 and <= 3, the L2 norm of the component is returned.
Examples
julia> using StaticArrays
julia> e = @SArray([1 4 7;2 5 8;3 6 9]);
julia> abs(e)
16.881943016134134
julia> abs(e,1)
3.7416573867739413
PlasmaEquilibriumToolkit.MagneticCoordinateGrid
— MethodMagneticCoordinateGrid(C::AbstractMagneticCoordinates,α::Union{AbstractArray,Real},β::Union{AbstractArray,Real},η::Union{AbstractArray,Real}
Generate a grid of magnetic coordiantes of type C
, returning a StructArray{C}
containing the coordinates
PlasmaEquilibriumToolkit.basis_vectors
— Methodbasis_vectors(b::BasisType,t::Transformation,x::AbstractArray{T},eq::AbstractMagneticEquilibrium) where T <: AbstractMagneticCoordinates
Generate basis vectors of type b
for the magnetic coordinates defined by the transformation t
at the points given by x
using data from the magnetic equilibrium eq
. The routine providing the point transformation designated by t
is defined in the respective equilibrium module.
Examples
Generate covariant basis vectors for points along a field line defined by PestCoordinates using a VMEC equilibrium
julia> using PlasmaEquilibriumToolkit, VMEC, NetCDF
julia> wout = NetCDF.open("wout.nc"); vmec, vmec_data = readVmecWout(wout);
julia> s = 0.5; vmec_s = VmecSurface(s,vmec);
julia> x = PestCoordinates(s*vmec.phi(1.0)/2π*vmec.signgs,0.0,-π:2π/10:π);
julia> v = VmecFromPest()(x,vmec_s);
julia> e_v = basis_vectors(Covariant(),CartesianFromVmec(),v,vmec_s)
10-element Array{StaticArrays.SArray{Tuple{3,3},Float64,2,9},1}:
⋮
PlasmaEquilibriumToolkit.curvature_components
— Methodcurvature_components(e::BasisVectors,gradB::CoordinateVector)
Computes the normal and geodesic curvature vectors for a stright field line coordinate system with basis vectors ∇X = e[:,1]
and ∇Y = e[:,2]
.
See also: normalCurvature
, geodesicCurvature
PlasmaEquilibriumToolkit.geodesic_curvature
— Methodgeodesic_curvature(b::CoordinateVector,gradB::CoordinateVector,gradX::CoordinateVector)
Computes the geodesic curvature component defined in straight fieldline coordinates (X,Y,Z) with basis vectors defined by (∇X,∇Y,∇Z) with the relation κ_g = -(B × ∇B) ⋅ ∇X/(B²|∇X|)
PlasmaEquilibriumToolkit.grad_B_projection
— Methodgrad_B_projection(e::BasisVectors,∇B::CoordinateVector)
Computes the projection of B × ∇B/B² onto the perpendicular coordinate vectors given by ∇X = e[:,1]
and ∇Y = e[:,2]
.
PlasmaEquilibriumToolkit.jacobian
— Methodjacobian(t::BasisType,e::BasisVectors{T})
Compute the Jacobian of the covariant/contravariatn basis vectors e
given by $J = e₁ ⋅ e₂ × e₃$ (covariant) or $J = 1/(e¹ ⋅ e² × e³)$ (contravariant).
PlasmaEquilibriumToolkit.metric
— Methodmetric2(e::BasisVectors)
Computes the metric tensor components gᵘᵛ
or gᵤᵥ
for a set of basis vectors e
in a 6 element SVector: [g11 = e[:,1]*e[:,1]
, g12 = e[:,1]*e[:,2]
,…]
PlasmaEquilibriumToolkit.normal_curvature
— MethodnormalCurvature(B::CoordinateVector,∇B::CoordinateVector,∇X::CoordinateVector,∇Y::CoordinateVector)
Computes the normal curvature component defined in straight fieldline coordinates (X,Y,Z) with basis vectors defined by (∇X,∇Y,∇Z) with the relation κₙ = (B × ∇B) ⋅ ((∇X⋅∇X)∇Y - (∇X⋅∇Y)∇X)/(B³|∇X|)
PlasmaEquilibriumToolkit.transform_basis
— Methodtransform_basis(t::Transformation,x::AbstracyArray{AbstractMagneticCoordinates},e::AbstractArray{BasisVectors{T}},eq::AbstractMagneticEquilibrium)
Perform a change of basis for magnetic coordinates denoted by the transformation t
, using the provided coordinate charts and basis vectors.
PlasmaEquilibriumToolkit.transform_basis
— Methodtransform_basis(t::BasisTransformation, e::BasisVectors)
transform_basis(t::BasisTransformation, e::BasisVectors, jacobian)
transform_basis(t::BasisTransformation, e::AbstractArray{BasisVectors})
transform_basis(t::BasisTransformation, e::AbstractArray{BasisVectors}, jacobian::AbstractArray)
Transform the covariant/contravariant basis given by e
to a contravariant/covariant basis by computing the transformation Jacobian from the basis vectors; for ContravariantFromCovariant()
the Jacobian is given by $J = e₁ ⋅ e₂ × e₃$ and for CovariantFromContravariant()
the Jacobian is given by $J = 1.0/(e¹ ⋅ e² × e³)$.