get_flux_c_ Function

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.

Type Bound

radbelt_type

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:

  • 1 -- particle species: electrons, solar activity: min
  • 2 -- particle species: electrons, solar activity: max
  • 3 -- particle species: protons, solar activity: min
  • 4 -- particle species: protons, solar activity: max

Return Value real(kind=wp)

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


Calls

proc~~get_flux_c_~~CallsGraph proc~get_flux_c_ radbelt_module::radbelt_type%get_flux_c_ proc~aep8 trmfun_module::trm_type%aep8 proc~get_flux_c_->proc~aep8 proc~igrfc shellig_module::shellig_type%igrfc proc~get_flux_c_->proc~igrfc proc~get_data_file_dir~2 trmfun_module::trm_type%get_data_file_dir proc~aep8->proc~get_data_file_dir~2 proc~trara1 trmfun_module::trm_type%trara1 proc~aep8->proc~trara1 proc~feldc shellig_module::shellig_type%feldc proc~igrfc->proc~feldc proc~feldcof shellig_module::shellig_type%feldcof proc~igrfc->proc~feldcof proc~findb0 shellig_module::shellig_type%findb0 proc~igrfc->proc~findb0 proc~shellc shellig_module::shellig_type%shellc proc~igrfc->proc~shellc proc~extrashc shellig_module::extrashc proc~feldcof->proc~extrashc proc~get_data_file_dir shellig_module::shellig_type%get_data_file_dir proc~feldcof->proc~get_data_file_dir proc~getshc shellig_module::getshc proc~feldcof->proc~getshc proc~intershc shellig_module::intershc proc~feldcof->proc~intershc proc~stoer shellig_module::shellig_type%stoer proc~findb0->proc~stoer proc~shellg shellig_module::shellig_type%shellg proc~shellc->proc~shellg proc~trara2 trmfun_module::trm_type%trara2 proc~trara1->proc~trara2 proc~shellg->proc~stoer proc~geo_to_cart shellig_module::geo_to_cart proc~shellg->proc~geo_to_cart proc~feldi shellig_module::shellig_type%feldi proc~stoer->proc~feldi

Called by

proc~~get_flux_c_~~CalledByGraph proc~get_flux_c_ radbelt_module::radbelt_type%get_flux_c_ none~get_flux radbelt_module::radbelt_type%get_flux none~get_flux->proc~get_flux_c_ proc~get_flux_c radbelt_module::get_flux_c proc~get_flux_c->none~get_flux proc~get_flux_g radbelt_module::get_flux_g proc~get_flux_g->none~get_flux proc~get_flux_g_c radbelt_c_module::get_flux_g_c proc~get_flux_g_c->none~get_flux interface~get_flux radbelt_module::get_flux interface~get_flux->proc~get_flux_c interface~get_flux->proc~get_flux_g

Source Code

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

        class(radbelt_type), intent(inout) :: me
        real(wp), dimension(3), intent(in) :: v
        real(wp), intent(in) :: year !! decimal year for which geomagnetic field is to
                                  !! be calculated (e.g.:1995.5 for day 185 of 1995)
        real(wp), intent(in) :: e !! minimum energy
        integer, intent(in) :: imname !! which method to use:
                                   !!
                                   !! * 1 -- particle species: electrons, solar activity: min
                                   !! * 2 -- particle species: electrons, solar activity: max
                                   !! * 3 -- particle species: protons, solar activity: min
                                   !! * 4 -- particle species: protons, solar activity: max
        real(wp) :: flux !! The flux of particles above the given energy, in units of cm^-2 s^-1.

        real(wp) :: xl !! l value
        real(wp) :: bbx

        call me%igrf%igrfc(v, year, xl, bbx)
        call me%trm%aep8(e, xl, bbx, imname, flux)

    end function get_flux_c_