le_grids_config_type Derived Type

type, public, extends(abstract_config_type) :: le_grids_config_type

Used to represent the input configuration of le_grids


Contents

Source Code


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 :: bouncefuzz = 10*epsilon(0.0)

Acts as a small tolerance in deciding if a particular pitch angle is forbidden at a particular theta grid point. Rather than considering particles to be forbidden if 1.0 - lambda * B(theta) < 0 we treat them as forbidden if 1.0 - lambda * B(theta) < -bouncefuzz. This provides a small "buffer" to account for floating round off in the calculation of lambda and B.

logical, public :: genquad = .false.

If true use generalised quadrature scheme for velocity integrals and energy grid. See G. Wilkie thesis.

integer, public :: negrid = -10

Sets the number of energy grid points to use. If specified then overrides values of nesuper = min(negrid/10 + 1, 4) and nesub = negrid - nesuper. If not set then the total number of energy grid points is just negrid = nesub + nesuper. The energy grid is split into two regions at a point controlled by vcut.

integer, public :: nesub = 8

Sets the number of energy grid points below the cutoff.

integer, public :: nesuper = 2

Sets the number of energy grid points above the cutoff.

logical, public :: new_trap_int = .false.

If true then use a more accurate integration method for the trapped pitch angle contribution to velocity integrals. The default uses a simple, but potentially less accurate, finite difference method.

logical, public :: new_trap_int_split = .false.

If true then split the trapped region into two symmetric regions when calculating the integration weights associated with the new_trap_int approach.

integer, public :: npassing = -1

The number of untrapped pitch-angles moving in one direction along field line.

integer, public :: ngauss = 5

The number of untrapped pitch-angles moving in one direction along field line is 2*ngauss if npassing is not set. Note that the number of trapped pitch angles is directly related to the number of theta grid points (per ) and is ntheta/2 + 1.

integer, public :: nmax = 500

Influences the minimum number of grid points in an integration subinterval used in calculating the integration grid weights.

logical, public :: radau_gauss_grid = .true.

The default lambda grid is now gauss-radau, for passing particles up to and including the passing-trapped boundary (wfb). The grid gives finite weight to the class of particles which bounce at theta = +/- pi (wfb) in the passing integral Note that the old gauss-legendre is used in the case that trapped_particles =.false.

logical, public :: split_passing_region = .false.

If true then we split the passing region into two separate regions for the purpose of choosing the integration weights. The split point is determined automatically by an attempt to minimise the passing weights error.

logical, public :: test = .false.

If true then just does a few tests and writes pitch angle and energy grids to screen before aborting the simulation.

logical, public :: trapped_particles = .true.

If set to false then the lambda grid weighting wl is set to zero for trapped pitch angles. This means that integrals over velocity space do not include a contribution from trapped particles which is equivalent to the situation where eps<=0.0. Trapped particle drifts are not set to zero so "trapped" particles still enter the source term through wdfac. At least for s-alpha the drifts are the main difference (please correct if not true) between the eps<=0.0 and the trapped_particles = .false. cases as Bmag is not a function of theta in the eps<=0.0 case whilst it is in the trapped_particles = .false. case.

real, public :: vcut = 2.5

No. of standard deviations from the standard Maxwellian beyond which the distribution function will be set to 0.

character(len=20), public :: wfbbc_option = 'default'

Set boundary condition for WFB in the linear/parallel solve:

  • "default", "mixed": The previous default boundary condition which mixes the passing and trapped boundary conditions
  • "passing": Treats the WFB using a passing particle boundary condition
  • "trapped": Treats the WFB using a trapped particle boundary condition
real, public :: wgt_fac = 10.0

Influences the maximum integration weight allowed.


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, 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_le_grids_config

  • private subroutine read_le_grids_config(self)

    Reads in the le_grids_knobs namelist and populates the member variables

    Arguments

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

procedure, public :: write => write_le_grids_config

  • private subroutine write_le_grids_config(self, unit)

    Writes out a namelist representing the current state of the config object

    Arguments

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

procedure, public :: reset => reset_le_grids_config

procedure, public :: broadcast => broadcast_le_grids_config

  • private subroutine broadcast_le_grids_config(self)

    Broadcasts all config parameters so object is populated identically on all processors

    Arguments

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

procedure, public, nopass :: get_default_name => get_default_name_le_grids_config

procedure, public, nopass :: get_default_requires_index => get_default_requires_index_le_grids_config

Source Code

  type, extends(abstract_config_type) :: le_grids_config_type
     ! namelist : le_grids_knobs
     ! indexed : false
     !> Acts as a small tolerance in deciding if a particular pitch
     !> angle is forbidden at a particular theta grid point. Rather
     !> than considering particles to be forbidden if `1.0 - lambda *
     !> B(theta) < 0` we treat them as forbidden if `1.0 - lambda *
     !> B(theta) < -bouncefuzz`. This provides a small "buffer" to
     !> account for floating round off in the calculation of `lambda`
     !> and `B`.
     real :: bouncefuzz = 10*epsilon(0.0)
     !> If true use generalised quadrature scheme for velocity
     !> integrals and energy grid. See [G. Wilkie
     !> thesis](https://drum.lib.umd.edu/handle/1903/17302).
     logical :: genquad = .false.
     !> Sets the number of energy grid points to use. If specified
     !> then overrides values of `nesuper = min(negrid/10 + 1,
     !> 4)` and `nesub = negrid - nesuper`. If not set then the total
     !> number of energy grid points is just `negrid = nesub +
     !> nesuper`. The energy grid is split into two regions at a point
     !> controlled by `vcut`.
     integer :: negrid = -10
     !> Sets the number of energy grid points below the cutoff.
     integer :: nesub = 8
     !> Sets the number of energy grid points above the cutoff.
     integer :: nesuper = 2
     !> If `true` then use a more accurate integration method for the
     !> trapped pitch angle contribution to velocity integrals. The
     !> default uses a simple, but potentially less accurate, finite
     !> difference method.
     logical :: new_trap_int = .false.
     !> If `true` then split the trapped region into two symmetric
     !> regions when calculating the integration weights associated
     !> with the new_trap_int approach.
     logical :: new_trap_int_split = .false.
     !> The number of untrapped pitch-angles moving in one direction
     !> along field line.
     integer :: npassing = -1
     !> The number of untrapped pitch-angles moving in one direction
     !> along field line is `2*ngauss` if [[le_grids_knobs:npassing]] is
     !> not set. Note that the number of trapped pitch angles is
     !> directly related to the number of theta grid points (per
     !> \(2\pi\)) and is `ntheta/2 + 1`.
     integer :: ngauss = 5
     !> Influences the minimum number of grid points in an integration
     !> subinterval used in calculating the integration grid weights.
     integer :: nmax = 500
     !> The default lambda grid is now gauss-radau, for passing
     !> particles up to and including the passing-trapped boundary
     !> (wfb). The grid gives finite weight to the class of particles
     !> which bounce at theta = +/- pi (wfb) in the passing integral
     !> Note that the old gauss-legendre is used in the case that
     !> trapped_particles =.false.
     !>
     !> @note In reference to above comment -
     !> code doesn't seem to change radau_gauss_grid when
     !> trapped_particles is false, should it?
     logical :: radau_gauss_grid = .true.
     !> If true then we split the passing region into two separate
     !> regions for the purpose of choosing the integration weights.
     !> The split point is determined automatically by an attempt
     !> to minimise the passing weights error.
     logical :: split_passing_region = .false.
     !> If `true` then just does a few tests and writes pitch angle
     !> and energy grids to screen before aborting the simulation.
     logical :: test = .false.
     !> If set to `false` then the lambda grid weighting `wl` is set
     !> to zero for trapped pitch angles. This means that integrals
     !> over velocity space do not include a contribution from trapped
     !> particles which is equivalent to the situation where
     !> `eps<=0.0`.  Trapped particle drifts are not set to zero so
     !> "trapped" particles still enter the source term through
     !> `wdfac`. At least for s-alpha the drifts are the main
     !> difference (*please correct if not true*) between the
     !> `eps<=0.0` and the `trapped_particles = .false.` cases as
     !> `Bmag` is not a function of theta in the `eps<=0.0` case
     !> whilst it is in the `trapped_particles = .false.` case.
     logical :: trapped_particles = .true.
     !> No. of standard deviations from the standard Maxwellian beyond
     !> which the distribution function will be set to 0.
     real :: vcut = 2.5
     !> Set boundary condition for WFB in the linear/parallel solve:
     !>
     !> - "default", "mixed": The previous default boundary condition which
     !>   mixes the passing and trapped boundary conditions
     !> - "passing": Treats the WFB using a passing particle boundary condition
     !> - "trapped": Treats the WFB using a trapped particle boundary condition
     character(len = 20) :: wfbbc_option = 'default'
     !> Influences the maximum integration weight allowed.
     real :: wgt_fac = 10.0
   contains
     procedure, public :: read => read_le_grids_config
     procedure, public :: write => write_le_grids_config
     procedure, public :: reset => reset_le_grids_config
     procedure, public :: broadcast => broadcast_le_grids_config
     procedure, public, nopass :: get_default_name => get_default_name_le_grids_config
     procedure, public, nopass :: get_default_requires_index => get_default_requires_index_le_grids_config
  end type le_grids_config_type