equilibrium_t Derived Type

type, public, abstract :: equilibrium_t

Class-instance variables


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(init), public, deferred, pass(self) :: init

  • subroutine init(self, filename) Prototype

    Initialises equilibrium, i.e.~reads any required parameters from file.

    Arguments

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

procedure(display), public, deferred, pass(self) :: display

  • subroutine display(self) Prototype

    Print to console minimal information about the equilibrium

    Arguments

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

procedure(debug), public, deferred, pass(self) :: debug

  • subroutine debug(self) Prototype

    Print to console extensive information about the equilibrium

    Arguments

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

procedure(is_axisymmetric), public, deferred, pass(self) :: is_axisymmetric

  • function is_axisymmetric(self) Prototype

    Returns if equilibrium is axisymmetric or 3D (stellarator)

    Arguments

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

    Return Value logical

procedure(rho), public, deferred, pass(self) :: rho

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

    Flux surface label normalised such that rho = 0 at the magnetic axis, and rho = 1 at the separatrix

    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(bx), public, deferred, pass(self) :: bx

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

    Radial field component normalized to on axis field strength.

    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(by), public, deferred, pass(self) :: by

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

    Vertical field component normalized to on axis field strength.

    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(btor), public, deferred, pass(self) :: btor

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

    Toroidal field component normalized to on axis field strength.

    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(jacobian), public, deferred, pass(self) :: jacobian

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

    Jacobian of underlying coordinate system

    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(epol), public, deferred, pass(self) :: epol

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

    Poloidal (along flux surface) unit vector.

    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)

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

procedure(erad), public, deferred, pass(self) :: erad

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

    Radial (across flux surface) unit vector.

    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)

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

procedure(district), public, deferred, pass(self) :: district

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

    Returns a descriptor which specifies whether a given point is in the closed, open or private flux region, and whether that point is within the computational domain

    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 integer

procedure(in_vessel), public, deferred, pass(self) :: in_vessel

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

    Returns a boolean of whether a given point is within the reactor vessel (inside the bounding surface of the first-wall and divertor) or outside this bounding surface. N.b. DISTRICT_CORE is still within the first-wall/divertor

    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 logical

procedure(mag_axis_loc), public, deferred, pass(self) :: mag_axis_loc

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

    Calculated the coordinates of the magnetic axis

    Arguments

    Type IntentOptional Attributes Name
    class(equilibrium_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