leq_type Derived Type

type, public, extends(abstract_geo_type) :: leq_type


Contents

Source Code


Components

Type Visibility Attributes Name Initial
character(len=:), public, allocatable :: type_name

Name of the specific implementation

logical, public :: initialised = .false.

Has this instance been initialised

logical, public :: has_full_theta_range = .false.

Does this implementation have the full or half theta grid? If only half then we assume up-down symmetry and use this to set gradients etc.

integer, public :: nr

Typically the number of radial and theta grid points.

integer, public :: nt

Typically the number of radial and theta grid points.

real, public :: psi_0
real, public :: psi_a
real, public :: B_T
real, public :: beta_0
real, public :: R_mag
real, public :: Z_mag
real, public :: aminor
real, public, allocatable, dimension(:,:) :: R_psi

2D map of co-ordinates to R, Z or B

real, public, allocatable, dimension(:,:) :: Z_psi

2D map of co-ordinates to R, Z or B

real, public, allocatable, dimension(:,:) :: B_psi

2D map of co-ordinates to R, Z or B

real, public, allocatable, dimension(:,:) :: eqth

2D maps of theta and eqpsi

real, public, allocatable, dimension(:,:) :: eqpsi_2d

2D maps of theta and eqpsi

real, public, allocatable, dimension (:,:,:) :: dpcart

Minor radius gradient in cylindrical coordinates

real, public, allocatable, dimension (:,:,:) :: dtcart

Theta gradient in cylindrical coordinates

real, public, allocatable, dimension (:,:,:) :: dbcart

B gradient in cylindrical coordinates

real, public, allocatable, dimension (:,:,:) :: dpbish

Minor radius gradient in Bishop coordinates

real, public, allocatable, dimension (:,:,:) :: dtbish

Theta gradient in Bishop coordinates

real, public, allocatable, dimension (:,:,:) :: dbbish

B gradient in Bishop coordinates

real, public, allocatable, dimension(:) :: eqpsi

Equilibrium psi (not normalised) and pressure grids

real, public, allocatable, dimension(:) :: pressure

Equilibrium psi (not normalised) and pressure grids

real, public, allocatable, dimension(:) :: psi_bar

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public, allocatable, dimension(:) :: fp

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public, allocatable, dimension(:) :: qsf

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public, allocatable, dimension(:) :: beta

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public, allocatable, dimension(:) :: diam

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public, allocatable, dimension(:) :: rc

Normalised psi, fp, q, beta, diameter and R_centre values on grid

real, public :: thetaShift

Used to translate definition of the poloidal angle such that the location of the maximium magnetic field is at pi (see section 3.2.3 of Ball MIT Masters thesis or section 2.4 of "GS2 analytic geometry specification")

type(flux_surface_type), public :: surf

Type-Bound Procedures

procedure, public, :: initialise

  • public subroutine initialise(self, inputs, psi_0_out, psi_a_out, rmaj, B_T0, avgrmid)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(inout) :: self
    type(geo_input_type), intent(in) :: inputs
    real, intent(out) :: psi_0_out
    real, intent(out) :: psi_a_out
    real, intent(out) :: rmaj
    real, intent(out) :: B_T0
    real, intent(out) :: avgrmid

procedure, public, :: gradient

Calculate the derivative of rp w.r.t. R and Z and return the modulus sqrt(drp/dR ^ 2 + drp/dZ^2). I.e. return |grad rp|. Parameters are: - rgrid: Submodule radial grid (r_circ for EFIT, rp otherwise) - theta: grid of theta - gradf(:,1): |grad psi| - gradf(:,2): 0.0 - char: if char = 'R', return |grad pressure| instead - char: if char = 'T', then return theta gradient in bishop form - gradf(:,1): (dtheta/dZ d rp/dR - dtheta/dR drp/dZ)/ |grad rp| - gradf(:,2): (dtheta/dR d rp/dR + dtheta/dZ drp/dZ)/ |grad rp| - rp: value of rp on the flux surface where we want the grad - nth: number of theta points - ntgrid: lower index of theta array is -ntgrid

  • public subroutine gradient(self, rgrid, theta, grad, char, rp, nth_used, ntm, use_bishop)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in), dimension (-ntm:) :: rgrid
    real, intent(in), dimension (-ntm:) :: theta
    real, intent(out), dimension (-ntm:,:) :: grad
    character(len=1), intent(in) :: char
    real, intent(in) :: rp
    integer, intent(in) :: nth_used
    integer, intent(in) :: ntm
    logical, intent(in) :: use_bishop

procedure, public, :: calculate_gradients

  • public pure subroutine calculate_gradients(self)

    Given theta, R, Z, B and psi on 2D grids calculate the index space derivatives in the two grid dimensions and use these to find gradients in cartesian and Bishop space.

    Arguments

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

procedure, public, :: rhofun => rhofun_null

  • public function rhofun_null(self, pbar) result(rhofun)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: derm

  • public pure subroutine derm(self, f, dfm, char)

    Calculate the derivative of f w.r.t to the radial and poloidal indexes (i.e. calculate the finite differences without dividing by delta psi and delta theta). - dfm(:,:,1) is the psi derivative - dfm(:,:,2) is the theta derivative - char specifies the periodicity in theta - 'E', 'O' mean continuous at theta = +/- pi - 'T' means a jump of 2 pi at theta = +/- pi

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in), dimension(:,:) :: f
    real, intent(out), dimension(:,:,:) :: dfm
    character(len=1), intent(in) :: char

procedure, public, :: invR

  • public function invR(self, r, theta)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta

    Return Value real

procedure, public, :: eqdbish

  • public pure subroutine eqdbish(self, dcart, dpcart, dbish)

    Convert gradients of a function f w.r.t. R,Z into bishop form. Note that dbish(1, :, :) is not valid

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in), dimension(:, :, :) :: dcart

    dcart(:,:,1) is gradient of f w.r.t. R; dcart(:,:,2) is gradient of f w.r.t. Z

    real, intent(in), dimension(:, :, :) :: dpcart

    dcart(:,:,1) is gradient of f w.r.t. R; dcart(:,:,2) is gradient of f w.r.t. Z

    real, intent(out), dimension(:, :, :) :: dbish

    dbish(:,:,1) is set to (df/dR dpsi/dR + df/dZ dpsi/dZ)/|grad psi|; dbish(:,:,2) is set to (-df/dR dpsi/dZ + df/dZ dpsi/dR)/|grad psi|. Note that in the special case where f=psi: dbish(:,:,1) is |grad psi|; dbish(:,:,2) is 0

procedure, public, :: eqdcart

  • public pure subroutine eqdcart(self, dfm, drm, dzm, dfcart)

    Converts derivatives w.r.t. (psi_index,theta_index) to derivatives w.r.t. (R,Z). Note that dfcart(1, :, :) is not valid

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in), dimension (:,:,:) :: dfm

    dfm(:,:,1) is deriv w.r.t. psi_index (i.e. (df/dpsi)_theta * delta psi); dfm(:,:,2) is deriv w.r.t. theta_index

    real, intent(in), dimension (:,:,:) :: drm

    dfm(:,:,1) is deriv w.r.t. psi_index (i.e. (df/dpsi)_theta * delta psi); dfm(:,:,2) is deriv w.r.t. theta_index

    real, intent(in), dimension (:,:,:) :: dzm

    dfm(:,:,1) is deriv w.r.t. psi_index (i.e. (df/dpsi)_theta * delta psi); dfm(:,:,2) is deriv w.r.t. theta_index

    real, intent(out), dimension (:,:,:) :: dfcart

    dfcart(:,:,1) is deriv w.r.t. R; dfcart(:,:,2) is deriv w.r.t. Z

procedure, public, :: dealloc_common_arrays

procedure, public, :: alloc_common_arrays

procedure, public, :: hahm_burrell

  • public subroutine hahm_burrell(self, a)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_geo_type), intent(in) :: self
    real, intent(in) :: a

procedure, public, :: dealloc_arrays

procedure, public, :: alloc_arrays

procedure, public, :: alloc_special_arrays => alloc_special_arrays_null

procedure, public, :: dealloc_special_arrays => dealloc_special_arrays_null

procedure, public, :: read_and_set

  • private subroutine read_and_set(self, inputs)

    Set module-level parameters and initialise the leq module

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(inout) :: self
    type(geo_input_type), intent(in) :: inputs

procedure, public, :: finalise

  • private pure subroutine finalise(self)

    Finalise leq module

    Arguments

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

procedure, public, :: finish_setup

  • private pure subroutine finish_setup(self)

    Finish the leq setup

    Arguments

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

procedure, public, :: qfun

  • private pure function qfun(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: dbtori

  • private pure function dbtori(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: btori

  • private pure function btori(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: diameter

  • private pure function diameter(self, rp)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: rp

    Return Value real

procedure, public, :: rcenter

  • private pure function rcenter(self, rp)

    Return the major radius of the centre of the flux surface

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: rp

    Return Value real

procedure, public, :: zpos => Zpos

  • private pure function Zpos(self, r, theta)

    Return the axial coordinate at

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta

    Return Value real

procedure, public, :: rpos => Rpos

  • private pure function Rpos(self, r, theta)

    Return the major radius at

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta

    Return Value real

procedure, public, :: eqitem

  • private pure function eqitem(self, r_in, theta_in, f, char) result(fstar)

    Calculates fstar which is f interpolated at the location theta_in

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r_in
    real, intent(in) :: theta_in

    Desired theta location

    real, intent(in), dimension(:, :) :: f

    Array to interpolate

    character(len=1), intent(in) :: char

    FIXME: Unused

    Return Value real

procedure, public, :: rfun

  • private pure function rfun(self, r, theta)

    FIXME : Add documentation | Should this be equivalent to psi as in others?

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta

    Return Value real

procedure, public, :: psi

  • private pure function psi(self, r, theta)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta

    Return Value real

procedure, public, :: betafun

  • private pure function betafun(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: dpfun

  • private pure function dpfun(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: pfun

  • private pure function pfun(self, pbar)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: pbar

    Return Value real

procedure, public, :: geo_FourierSeries

  • private pure subroutine geo_FourierSeries(self, r, theta, Rpos, Zpos)

    Calculate the cylindrical at for the Fourier model

    See section 5.1.1 of Ball Oxford PhD thesis or section 2.1.4 of "GS2 analytic geometry specification"

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta
    real, intent(out) :: Rpos
    real, intent(out) :: Zpos

procedure, public, :: geo_generalizedElongation

  • private pure subroutine geo_generalizedElongation(self, r, theta, Rpos, Zpos)

    Calculate the cylindrical at for the generalised elongation model

    See section 5.1.2 of Ball Oxford PhD thesis or section 2.1.3 of "GS2 analytic geometry specification"

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta
    real, intent(out) :: Rpos
    real, intent(out) :: Zpos

procedure, public, :: geo_global

  • private pure subroutine geo_global(self, r, theta, Rpos, Zpos)

    Calculate the cylindrical at for the global MHD model

    See section 2.1.2 of "GS2 analytic geometry specification"

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta
    real, intent(out) :: Rpos
    real, intent(out) :: Zpos

procedure, public, :: geo_Miller

  • private pure subroutine geo_Miller(self, r, theta, Rpos, Zpos)

    Calculate the cylindrical (r, \theta) for the generalised Miller model

    See section 3.2.1 of Ball MIT Masters thesis or section 2.1.1 of "GS2 analytic geometry specification"

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta
    real, intent(out) :: Rpos
    real, intent(out) :: Zpos

procedure, public, :: geo_MillerExtendedHarmonic

  • private pure subroutine geo_MillerExtendedHarmonic(self, r, theta, Rpos, Zpos)

    Calculate the cylindrical at for the Miller Extended Harmonic model

    See PPCF 63 (2021) 012001 (5pp)

    Arguments

    Type IntentOptional Attributes Name
    class(leq_type), intent(in) :: self
    real, intent(in) :: r
    real, intent(in) :: theta
    real, intent(out) :: Rpos
    real, intent(out) :: Zpos

Source Code

  type, extends(abstract_geo_type) :: leq_type
     !> Used to translate definition of the poloidal angle such that the location
     !> of the maximium magnetic field is at pi (see section 3.2.3 of Ball MIT
     !> Masters thesis or section 2.4 of "GS2 analytic geometry specification")
     real :: thetaShift
     type(flux_surface_type) :: surf
   contains
     procedure :: finish_setup, finalise, read_and_set
     procedure :: rpos, zpos, rcenter, diameter, btori, dbtori, qfun
     procedure :: pfun, dpfun, betafun, psi, rfun, eqitem
     procedure :: geo_Miller, geo_global, geo_generalizedElongation, geo_FourierSeries
     procedure :: geo_MillerExtendedHarmonic
  end type leq_type