salpha_t Derived Type

type, public, extends(equilibrium_t) :: salpha_t

Type implementing an salpha like equilibrium. The equilibrium has concentric flux surfaces in a toroidal geometry. Shear as well as the q-factor can be specified. To get the correct poloidal flux in Wb = V*s the reference magnetic field strength has to be provided in Tesla and the reference length scale has to be provided in meters.


Components

Type Visibility Attributes Name Initial
logical, public :: initialized = .false.
real(kind=FP), public :: x0

Magnetic axis x = R/R0 (in normalised units)

real(kind=FP), public :: y0

Magnetic axis y = Z/R0 (in normalised units)

real(kind=FP), public :: phi0 = 0.0_FP

Magnetic axis phi

real(kind=FP), public :: xmin

Box limits

real(kind=FP), public :: xmax

Box limits

real(kind=FP), public :: ymin

Box limits

real(kind=FP), public :: ymax

Box limits

real(kind=FP), public :: rhomax

Global limits for rho (rho = normalised psi, n.b. there may also be region-specific limits defined in equi)

real(kind=FP), public :: rhomin

Global limits for rho (rho = normalised psi, n.b. there may also be region-specific limits defined in equi)

real(kind=FP), public :: rho_limiter

Radial position (start) of limiter


Type-Bound Procedures

procedure, public, pass(self) :: absb

  • private function absb(self, x, y, phi)

    Absolute value of magnetic field.

    Arguments

    Type IntentOptional Attributes Name
    class(equilibrium_t), intent(in) :: self
    real(kind=FP), intent(in) :: x

    3D position (x and y normalized)

    real(kind=FP), intent(in) :: y

    3D position (x and y normalized)

    real(kind=FP), intent(in) :: phi

    3D position (x and y normalized)

    Return Value real(kind=fp)

procedure, public, pass(self) :: bpol

  • private function bpol(self, x, y, phi)

    Magnetic field component b poloidal normalised to absolute value of B (on axis)

    Arguments

    Type IntentOptional Attributes Name
    class(equilibrium_t), intent(in) :: self
    real(kind=FP), intent(in) :: x

    3D position (x and y normalized)

    real(kind=FP), intent(in) :: y

    3D position (x and y normalized)

    real(kind=FP), intent(in) :: phi

    3D position (x and y normalized)

    Return Value real(kind=fp)

procedure, public :: init

  • private subroutine init(self, filename)

    Reads input parameters and initializes the equilibrium

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(inout) :: self
    character(len=*), intent(in), optional :: filename

    Filename, where parameters are read from

procedure, public :: display

  • private subroutine display(self)

    Prints to console information about the equilibrium

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

procedure, public :: debug

  • private subroutine debug(self)

    Prints to console information about the equilibrium

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

procedure, public :: is_axisymmetric

  • private function is_axisymmetric(self)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

    Return Value logical

procedure, public :: rho

  • private function rho(self, x, y, phi)

    Returns the flux surface label

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value real(kind=fp)

procedure, public :: psi

  • private function psi(self, x, y)

    Returns the poloidal flux function in SI units

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y

    Return Value real(kind=fp)

procedure, public :: bx

  • private function bx(self, x, y, phi)

    Returns magnetic field component bx normalised to B_0 (on axis)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value real(kind=fp)

procedure, public :: by

  • private function by(self, x, y, phi)

    Return the magnetic field component by normalised to B_0 (on axis)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

    Instnace of the type

    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value real(kind=fp)

procedure, public :: btor

  • private function btor(self, x, y, phi)

    Returns the magnetic field component btor normalised to B_0 (on axis)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value real(kind=fp)

procedure, public :: jacobian

  • private function jacobian(self, x, y, phi)

    Returns the jacobian of the geometry. The geometry is cylindrical and thus J = x

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value real(kind=fp)

procedure, public :: epol

  • private subroutine epol(self, x, y, phi, epolx, epoly)

    Calculates the unit vector along poloidal (along flux surface) direction

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi
    real(kind=FP), intent(out) :: epolx

    Unit vector component in x-direction

    real(kind=FP), intent(out) :: epoly

    Unit vector component in y-direction

procedure, public :: erad

  • private subroutine erad(self, x, y, phi, eradx, erady)

    Calculates unit vector along radial (across flux surface) direction

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi
    real(kind=FP), intent(out) :: eradx

    Unit vector component in x-direction

    real(kind=FP), intent(out) :: erady

    Unit vector component in y-direction

procedure, public :: district

  • private function district(self, x, y, phi)

    Returns in which district point (x, y) is (see module descriptors_m)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y
    real(kind=FP), intent(in) :: phi

    Return Value integer

procedure, public :: in_vessel

  • private function in_vessel(self, x, y, phi)

    Returns a boolean of whether a given point is within theta1 and theta2 (angular extent of limiter)

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: x

    x-coordinate

    real(kind=FP), intent(in) :: y

    y-coordinate

    real(kind=FP), intent(in) :: phi

    toroidal coordinate (unused)

    Return Value logical

procedure, public :: mag_axis_loc

  • private subroutine mag_axis_loc(self, phi, axis_x, axis_y)

    Returns the coordinates of magnetic axis

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

    Instance of class

    real(kind=FP), intent(in) :: phi

    Toroidal angle

    real(kind=FP), intent(out) :: axis_x

    x-coordinate of the magnetic axis

    real(kind=FP), intent(out) :: axis_y

    y-coordinate of the magnetic axis

procedure, public :: qval

  • private function qval(self, rho_pt)

    Returns the safety factor

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self
    real(kind=FP), intent(in) :: rho_pt

    Flux surface label

    Return Value real(kind=fp)

procedure, public, nopass :: theta

  • private function theta(x, y)

    Returns the geometric poloidal angle

    Arguments

    Type IntentOptional Attributes Name
    real(kind=FP), intent(in) :: x
    real(kind=FP), intent(in) :: y

    Return Value real(kind=fp)

procedure, public :: x_point_psi

  • private function x_point_psi(self)

    Return the value of psi at last closed flux surface

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

    Return Value real(kind=fp)

procedure, public :: o_point_psi

  • private function o_point_psi(self)

    Return the value of psi at magnetic axis

    Arguments

    Type IntentOptional Attributes Name
    class(salpha_t), intent(in) :: self

    Return Value real(kind=fp)