radbelt_module Module

Main module.

See also

  • https://ccmc.gsfc.nasa.gov/pub/modelweb/geomagnetic/igrf/fortran_code/bilcal.for
  • https://ccmc.gsfc.nasa.gov/pub/modelweb/radiation_belt/radbelt/fortran_code/radbelt.for

Uses

  • module~~radbelt_module~~UsesGraph module~radbelt_module radbelt_module module~radbelt_kinds_module radbelt_kinds_module module~radbelt_module->module~radbelt_kinds_module module~shellig_module shellig_module module~radbelt_module->module~shellig_module module~trmfun_module trmfun_module module~radbelt_module->module~trmfun_module iso_fortran_env iso_fortran_env module~radbelt_kinds_module->iso_fortran_env module~shellig_module->module~radbelt_kinds_module module~trmfun_module->module~radbelt_kinds_module

Used by

  • module~~radbelt_module~~UsedByGraph module~radbelt_module radbelt_module module~radbelt_c_module radbelt_c_module module~radbelt_c_module->module~radbelt_module

Interfaces

public interface get_flux

simple function versions for testing

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

    Calculate the flux of trapped particles at a specific location and time. This is just a function version of the class method from radbelt_type.

    Note

    This routine is not efficient at all since it will reload all the files every time it is called.

    Arguments

    Type IntentOptional Attributes Name
    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(v, year, e, imname) result(flux)

    Calculate the flux of trapped particles at a specific location and time. This is just a function version of the class method from radbelt_type.

    Note

    This routine is not efficient at all since it will reload all the files every time it is called.

    Arguments

    Type IntentOptional Attributes Name
    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.


Derived Types

type, public ::  radbelt_type

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

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_
procedure, private :: get_flux_c_
procedure, private :: get_flux_g_
procedure, public :: set_data_files_paths
procedure, public :: set_igrf_file_path
procedure, public :: set_trm_file_path

Functions

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_g(lon, lat, height, year, e, imname) result(flux)

Calculate the flux of trapped particles at a specific location and time. This is just a function version of the class method from radbelt_type.

Read more…

Arguments

Type IntentOptional Attributes Name
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.

public function get_flux_c(v, year, e, imname) result(flux)

Calculate the flux of trapped particles at a specific location and time. This is just a function version of the class method from radbelt_type.

Read more…

Arguments

Type IntentOptional Attributes Name
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.


Subroutines

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

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

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