theta_grid_gridgen_config_type Derived Type

type, public, extends(abstract_config_type) :: theta_grid_gridgen_config_type

Used to represent the input configuration of theta_grid_gridgen


Contents


Components

Type Visibility Attributes Name Initial
logical, public :: exist = .false.

Does the related namelist exist in the target input file?

integer, public :: index = 0

Used to hold the specific index of numbered namelists

logical, public :: skip_read = .false.

Do we want to skip the read step in init?

logical, public :: skip_broadcast = .false.

Do we want to skip the broadcast step in init?

real, public :: alknob = 0.0

Relative weighting of pitch-angle metric to metric

real, public :: bpknob = 1.e-8

Consider when the right grid point is equal to the target bmag.

FIXME: What does this mean?

real, public :: deltaw = 0.0

Parameter for weighted resolution in theta. Each theta grid point contributes to the resolution metric according to the function which has peaks at theta = +/- thetamax and deltaw controls the relative weighting of the theta dependent contribution.

real, public :: epsknob = 1e-5

Maximum difference between neighbouring points for determining if a point is an extremum.

real, public :: extrknob = 0.0

Used to set a "bonus" contribtion to resolution at B extrema with an even number of theta grid points. Those with an odd number of points and the assumed extrema at -pi have a metric of 1e20. Here extrknob can be used to bias the algorithm towards keeping extrema with an even number of points.

integer, public :: npadd = 2

Number of points between original grid points to oversample by

real, public :: regrid_tension = -1.0

Tension to use in interpolating splines for regrid of geometrical quantities Defaults to tension if not set.

logical, public :: skip_gridgen = .false.

If true then skip gridgen call and instead just use the existing grid.

real, public :: tension = 1.0

Tension for spline

real, public :: thetamax = 0.0

Parameter for weighted resolution in theta. Each theta grid point contributes to the resolution metric according to the function which has peaks at theta = +/- thetamax.

real, public :: widthw = 1.0

Parameter for weighted resolution in theta. Each theta grid point contributes to the resolution metric according to the function which has peaks at theta = +/- thetamax and widthw controls the scale length of the peaks.


Type-Bound Procedures

procedure, public, :: is_initialised => is_initialised_generic

procedure, public, :: init => init_generic

  • private subroutine init_generic(self, name, requires_index, index)

    Fully initialise the config object

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(inout) :: self
    character(len=*), intent(in), optional :: name
    logical, intent(in), optional :: requires_index
    integer, intent(in), optional :: index

procedure, public, :: setup => setup_generic

  • private subroutine setup_generic(self, name, requires_index, index)

    Do some standard setup/checking

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(inout) :: self
    character(len=*), intent(in), optional :: name
    logical, intent(in), optional :: requires_index
    integer, intent(in), optional :: index

procedure, public, :: write_namelist_header

  • private subroutine write_namelist_header(self, unit)

    Write the namelist header for this instance

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(in) :: self
    integer, intent(in) :: unit

procedure, public, :: get_name => get_name_generic

  • private function get_name_generic(self)

    Returns the namelist name. Not very useful at the moment but may want to do more interesting things in the future

    Arguments

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

    Return Value character(len=CONFIG_MAX_NAME_LEN)

procedure, public, :: get_requires_index => get_requires_index_generic

  • private function get_requires_index_generic(self)

    Returns the requires_index value. Allows access whilst keeping the variable private

    Arguments

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

    Return Value logical

procedure, public, :: set_smart_defaults => set_smart_defaults_null

  • private subroutine set_smart_defaults_null(self)

    An no-op implementation of the set_smart_defaults method. Unless over-ridden the specific config instance will have no smart defaults applied.

    Arguments

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

    Has to be intent in out as over-riding procedures need to change self

procedure, public, nopass :: write_namelist_footer

  • private subroutine write_namelist_footer(unit)

    Write the namelist footer

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: unit
  • private subroutine write_key_val_string(key, val, unit)

    Writes a {key,val} pair where the value is of type character

    Arguments

    Type IntentOptional Attributes Name
    character(len=*), intent(in) :: key
    character(len=*), intent(in) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_real(key, val, unit)

    Writes a {key,val} pair where the value is of type real

    Arguments

    Type IntentOptional Attributes Name
    character(len=*), intent(in) :: key
    real, intent(in) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_complex(key, val, unit)

    Writes a {key,val} pair where the value is of type complex

    Arguments

    Type IntentOptional Attributes Name
    character(len=*), intent(in) :: key
    complex, intent(in) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_integer(key, val, unit)

    Writes a {key,val} pair where the value is of type integer

    Arguments

    Type IntentOptional Attributes Name
    character(len=*), intent(in) :: key
    integer, intent(in) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_logical(key, val, unit)

    Writes a {key,val} pair where the value is of type logical

    Arguments

    Type IntentOptional Attributes Name
    character(len=*), intent(in) :: key
    logical, intent(in) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_real_array(self, key, val, unit)

    Writes a {key,val} pair where the value is of type real array

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(in) :: self
    character(len=*), intent(in) :: key
    real, intent(in), dimension(:) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_complex_array(self, key, val, unit)

    Writes a {key,val} pair where the value is of type complex array

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(in) :: self
    character(len=*), intent(in) :: key
    complex, intent(in), dimension(:) :: val
    integer, intent(in) :: unit
  • private subroutine write_key_val_integer_array(self, key, val, unit)

    Writes a {key,val} pair where the value is of type integer array

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_config_type), intent(in) :: self
    character(len=*), intent(in) :: key
    integer, intent(in), dimension(:) :: val
    integer, intent(in) :: unit

procedure, public :: read => read_theta_grid_gridgen_config

procedure, public :: write => write_theta_grid_gridgen_config

procedure, public :: reset => reset_theta_grid_gridgen_config

procedure, public :: broadcast => broadcast_theta_grid_gridgen_config

procedure, public, nopass :: get_default_name => get_default_name_theta_grid_gridgen_config

procedure, public, nopass :: get_default_requires_index => get_default_requires_index_theta_grid_gridgen_config

Source Code

  type, extends(abstract_config_type) :: theta_grid_gridgen_config_type
     ! namelist : theta_grid_gridgen_knobs
     ! indexed : false
     !> Relative weighting of pitch-angle metric to \(\theta\) metric
     real :: alknob = 0.0
     !> Consider when the right grid point is equal to the target bmag.
     !>
     !> FIXME: What does this mean?
     real :: bpknob = 1.e-8
     !> Parameter for weighted resolution in theta. Each theta grid point
     !> contributes to the resolution metric according to the function
     !> $$ 1 + deltaw *[ 1 / (1 + [theta-thetamax]**2/widthw**2) +
     !>                  1 / (1 + [theta+thetamax]**2/widthw**2)] $$
     !> which has peaks at theta = +/- thetamax and deltaw controls the
     !> relative weighting of the theta dependent contribution.
     real :: deltaw = 0.0
     !> Maximum difference between neighbouring points for determining
     !> if a point is an extremum.
     real :: epsknob = 1e-5
     !> Used to set a "bonus" contribtion to resolution at B extrema with an
     !> even number of theta grid points. Those with an odd number of points
     !> and the assumed extrema at -pi have a metric of 1e20. Here extrknob
     !> can be used to bias the algorithm towards keeping extrema with an
     !> even number of points.
     real :: extrknob = 0.0
     !> Number of points between original grid points to oversample \(\theta\) by
     integer :: npadd = 2
     !> Tension to use in interpolating splines for regrid of geometrical quantities
     !> Defaults to [[theta_grid_gridgen_config_type:tension]] if not set.
     real :: regrid_tension = -1.0
     !> If true then skip gridgen call and instead just use the existing grid.
     !>
     !> @note This is primarily for debugging and testing. When active we are
     !> effectively forced to assumed that the input bmag is symmetric and
     !> monotonic. If this is not true then we should not trust the results.
     logical :: skip_gridgen = .false.
     !> Tension for \(B\) spline
     real :: tension = 1.0
     !> Parameter for weighted resolution in theta. Each theta grid point
     !> contributes to the resolution metric according to the function
     !> $$ 1 + deltaw *[ 1 / (1 + [theta-thetamax]**2/widthw**2) +
     !>                  1 / (1 + [theta+thetamax]**2/widthw**2)] $$
     !> which has peaks at theta = +/- thetamax.
     real :: thetamax = 0.0
     !> Parameter for weighted resolution in theta. Each theta grid point
     !> contributes to the resolution metric according to the function
     !> $$ 1 + deltaw *[ 1 / (1 + [theta-thetamax]**2/widthw**2) +
     !>                  1 / (1 + [theta+thetamax]**2/widthw**2)] $$
     !> which has peaks at theta = +/- thetamax and widthw controls the
     !> scale length of the peaks.
     real :: widthw = 1.0
   contains
     procedure, public :: read => read_theta_grid_gridgen_config
     procedure, public :: write => write_theta_grid_gridgen_config
     procedure, public :: reset => reset_theta_grid_gridgen_config
     procedure, public :: broadcast => broadcast_theta_grid_gridgen_config
     procedure, public, nopass :: get_default_name => get_default_name_theta_grid_gridgen_config
     procedure, public, nopass :: get_default_requires_index => get_default_requires_index_theta_grid_gridgen_config
  end type theta_grid_gridgen_config_type