circular_toroidal_t Derived Type

type, public, extends(equilibrium_t) :: circular_toroidal_t

Type implementing a consistent analytical toroidal equilibrium with nested circular fluxsurfaces


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)


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)

    Initialises the toroidal equilibrium

    Arguments

    Type IntentOptional Attributes Name
    class(circular_toroidal_t), intent(inout) :: self

    Instance of class

    character(len=*), intent(in), optional :: filename

procedure, public :: display

  • private subroutine display(self)

    Prints to console information about the equilibrium and the boundary

    Arguments

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

    Instance of class

procedure, public :: debug

  • private subroutine debug(self)

    Prints to console extended debug information about the equilibrium

    Arguments

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

    Instance of class

procedure, public :: is_axisymmetric

  • private function is_axisymmetric(self)

    Returns if the equilibrium is axisymmetric or not

    Arguments

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

    Instance of class

    Return Value logical

procedure, public :: rho

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

    Returns the radial coordinate, i.e. distance to torus axis

    Arguments

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

    Instance of class

    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 :: bx

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

    Returns the x-component of the magnetic field

    Arguments

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

    Instance of class

    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 :: by

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

    Returns the y-component of the magnetic field

    Arguments

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

    Instance of class

    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 :: btor

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

    Returns the toroidal component btor of the magnetic field

    Arguments

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

    Instance of class

    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 :: jacobian

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

    Returns the Jacobian of the geometry

    Arguments

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

    Instance of class

    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 :: epol

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

    Calculates the the x- and y-component of the poloidal unit vector epol (along flux surfaces)

    Arguments

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

    Instance of class

    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)

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

    Components of poloidal unit vector erad

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

    Components of poloidal unit vector erad

procedure, public :: erad

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

    Calculates the the x- and y-component of the radial unit vector erad (across flux surfaces)

    Arguments

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

    Instance of class

    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)

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

    Components of radial unit vector erad

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

    Components of radial unit vector erad

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(circular_toroidal_t), intent(in) :: self

    Instance of class

    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 integer

procedure, public :: in_vessel

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

    Returns whether a given point is within the vessel or not

    Arguments

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

    Instance of class

    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 logical

procedure, public :: mag_axis_loc

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

    Calculated the coordinates of the magnetic axis

    Arguments

    Type IntentOptional Attributes Name
    class(circular_toroidal_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(circular_toroidal_t), intent(in) :: self

    Instance of class

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

    rho-coordinate (along minor radius)

    Return Value real(kind=fp)

procedure, public :: theta

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

    Returns the polar angle in polar coordinates

    Arguments

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

    Instance of class

    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)