get_flux_g_c Subroutine

public subroutine get_flux_g_c(ipointer, lon, lat, height, year, e, imname, flux) bind(C, name="get_flux_g_c")

C interface to get_flux_g.

Arguments

Type IntentOptional Attributes Name
integer(kind=c_intptr_t), intent(in) :: ipointer
real(kind=c_double), intent(in) :: lon

geodetic longitude in degrees (east)

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

geodetic latitude in degrees (north)

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

altitude in km above sea level

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

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

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

minimum energy

integer(kind=c_int), 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(kind=c_double), intent(out) :: flux

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


Calls

proc~~get_flux_g_c~~CallsGraph proc~get_flux_g_c radbelt_c_module::get_flux_g_c none~get_flux radbelt_module::radbelt_type%get_flux proc~get_flux_g_c->none~get_flux proc~int_pointer_to_f_pointer radbelt_c_module::int_pointer_to_f_pointer proc~get_flux_g_c->proc~int_pointer_to_f_pointer proc~get_flux_c_ radbelt_module::radbelt_type%get_flux_c_ none~get_flux->proc~get_flux_c_ proc~get_flux_g_ radbelt_module::radbelt_type%get_flux_g_ none~get_flux->proc~get_flux_g_ 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_flux_g_->proc~aep8 proc~igrf shellig_module::shellig_type%igrf proc~get_flux_g_->proc~igrf 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~feldcof shellig_module::shellig_type%feldcof proc~igrf->proc~feldcof proc~feldg shellig_module::shellig_type%feldg proc~igrf->proc~feldg proc~findb0 shellig_module::shellig_type%findb0 proc~igrf->proc~findb0 proc~shellg shellig_module::shellig_type%shellg proc~igrf->proc~shellg proc~feldc shellig_module::shellig_type%feldc proc~igrfc->proc~feldc proc~igrfc->proc~feldcof 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~shellc->proc~shellg proc~geo_to_cart shellig_module::geo_to_cart proc~shellg->proc~geo_to_cart proc~shellg->proc~stoer proc~trara2 trmfun_module::trm_type%trara2 proc~trara1->proc~trara2 proc~feldi shellig_module::shellig_type%feldi proc~stoer->proc~feldi

Source Code

    subroutine get_flux_g_c(ipointer, lon, lat, height, year, e, imname, flux) bind(C, name="get_flux_g_c")

        integer(c_intptr_t), intent(in) :: ipointer
        real(c_double), intent(in) :: lon !! geodetic longitude in degrees (east)
        real(c_double), intent(in) :: lat !! geodetic latitude in degrees (north)
        real(c_double), intent(in) :: height !! altitude in km above sea level
        real(c_double), intent(in) :: year !! decimal year for which geomagnetic field is to
                                      !! be calculated (e.g.:1995.5 for day 185 of 1995)
        real(c_double), intent(in) :: e !! minimum energy
        integer(c_int), 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(c_double), intent(out) :: flux !! The flux of particles above the given energy, in units of cm^-2 s^-1.

        type(radbelt_type), pointer :: p

        call int_pointer_to_f_pointer(ipointer, p)

        if (associated(p)) then

            flux = p%get_flux(lon, lat, height, year, e, imname)

        else
            error stop 'error in get_flux_g_c: class is not associated'
        end if

    end subroutine get_flux_g_c