radbelt_type Derived Type

type, public :: radbelt_type

the main class that can be used to get the flux.


Inherits

type~~radbelt_type~~InheritsGraph type~radbelt_type radbelt_type type~shellig_type shellig_type type~radbelt_type->type~shellig_type igrf type~trm_type trm_type type~radbelt_type->type~trm_type trm

Components

Type Visibility Attributes Name Initial
type(trm_type), private :: trm
type(shellig_type), private :: igrf

Type-Bound Procedures

generic, public :: get_flux => get_flux_g_, get_flux_c_

  • public function get_flux_g_(me, lon, lat, height, year, e, imname) result(flux)

    Calculate the flux of trapped particles at a specific location and time.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    real(kind=wp), intent(in) :: lon

    geodetic longitude in degrees (east)

    real(kind=wp), intent(in) :: lat

    geodetic latitude in degrees (north)

    real(kind=wp), intent(in) :: height

    altitude in km above sea level

    real(kind=wp), intent(in) :: year

    decimal year for which geomagnetic field is to be calculated (e.g.:1995.5 for day 185 of 1995)

    real(kind=wp), intent(in) :: e

    minimum energy

    integer, intent(in) :: imname

    which method to use:

    Read more…

    Return Value real(kind=wp)

    The flux of particles above the given energy, in units of cm^-2 s^-1.

  • public function get_flux_c_(me, v, year, e, imname) result(flux)

    Calculate the flux of trapped particles at a specific location and time. This is an alternate version of get_flux_g_ for cartesian coordinates.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    real(kind=wp), intent(in), dimension(3) :: v
    real(kind=wp), intent(in) :: year

    decimal year for which geomagnetic field is to be calculated (e.g.:1995.5 for day 185 of 1995)

    real(kind=wp), intent(in) :: e

    minimum energy

    integer, intent(in) :: imname

    which method to use:

    Read more…

    Return Value real(kind=wp)

    The flux of particles above the given energy, in units of cm^-2 s^-1.

procedure, private :: get_flux_c_

  • public function get_flux_c_(me, v, year, e, imname) result(flux)

    Calculate the flux of trapped particles at a specific location and time. This is an alternate version of get_flux_g_ for cartesian coordinates.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    real(kind=wp), intent(in), dimension(3) :: v
    real(kind=wp), intent(in) :: year

    decimal year for which geomagnetic field is to be calculated (e.g.:1995.5 for day 185 of 1995)

    real(kind=wp), intent(in) :: e

    minimum energy

    integer, intent(in) :: imname

    which method to use:

    Read more…

    Return Value real(kind=wp)

    The flux of particles above the given energy, in units of cm^-2 s^-1.

procedure, private :: get_flux_g_

  • public function get_flux_g_(me, lon, lat, height, year, e, imname) result(flux)

    Calculate the flux of trapped particles at a specific location and time.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    real(kind=wp), intent(in) :: lon

    geodetic longitude in degrees (east)

    real(kind=wp), intent(in) :: lat

    geodetic latitude in degrees (north)

    real(kind=wp), intent(in) :: height

    altitude in km above sea level

    real(kind=wp), intent(in) :: year

    decimal year for which geomagnetic field is to be calculated (e.g.:1995.5 for day 185 of 1995)

    real(kind=wp), intent(in) :: e

    minimum energy

    integer, intent(in) :: imname

    which method to use:

    Read more…

    Return Value real(kind=wp)

    The flux of particles above the given energy, in units of cm^-2 s^-1.

procedure, public :: set_data_files_paths

  • public subroutine set_data_files_paths(me, aep8_dir, igrf_dir)

    Set the paths to the data files. If not used or blank, the folder data/aep8 and data/igrf in the current working directory is assumed

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    character(len=*), intent(in) :: aep8_dir
    character(len=*), intent(in) :: igrf_dir

procedure, public :: set_igrf_file_path

  • public subroutine set_igrf_file_path(me, dir)

    Set the igrf path.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    character(len=*), intent(in) :: dir

procedure, public :: set_trm_file_path

  • public subroutine set_trm_file_path(me, dir)

    Set the trm path.

    Arguments

    Type IntentOptional Attributes Name
    class(radbelt_type), intent(inout) :: me
    character(len=*), intent(in) :: dir

Source Code

    type, public :: radbelt_type
    !! the main class that can be used to get the flux.
        private
        type(trm_type) :: trm
        type(shellig_type) :: igrf
    contains
        private
        generic, public :: get_flux => get_flux_g_, get_flux_c_
        procedure :: get_flux_g_, get_flux_c_
        procedure, public :: set_trm_file_path, &
            set_igrf_file_path, &
            set_data_files_paths
    end type radbelt_type