igrfc Subroutine

private subroutine igrfc(me, v, year, xl, bbx)

Alternate version of igrf for cartesian coordinates.

Type Bound

shellig_type

Arguments

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

cartesian coordinates in earth radii (6371.2 km) x-axis pointing to equator at 0 longitude y-axis pointing to equator at 90 long. z-axis pointing to north pole

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(out) :: xl

l-value

real(kind=wp), intent(out) :: bbx

b_total / b_equatorial ratio


Calls

proc~~igrfc~~CallsGraph proc~igrfc shellig_module::shellig_type%igrfc 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~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~~igrfc~~CalledByGraph proc~igrfc shellig_module::shellig_type%igrfc proc~get_flux_c_ radbelt_module::radbelt_type%get_flux_c_ proc~get_flux_c_->proc~igrfc 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

    subroutine igrfc(me, v, year, xl, bbx)

        class(shellig_type), intent(inout) :: me
        real(wp), dimension(3), intent(in) :: v !! cartesian coordinates in earth radii (6371.2 km)
                                             !! x-axis pointing to equator at 0 longitude
                                             !! y-axis pointing to equator at 90 long.
                                             !! z-axis pointing to north pole
        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(out) :: xl !! l-value
        real(wp), intent(out) :: bbx !! b_total / b_equatorial ratio

        real(wp) :: bab1, bdel, beq, bequ, dimo, rr0
        integer :: icode
        logical :: val
        real(wp), dimension(3) :: b

        real(wp), parameter :: stps = 0.05_wp

        ! JW : do we need to reset some or all of these ?
        me%sp = 0.0_wp
        me%xi = 0.0_wp
        me%h = 0.0_wp
        me%step = 0.20_wp
        me%steq = 0.03_wp

        call me%feldcof(year, dimo)
        call me%feldc(v, b)
        call me%shellc(v, dimo, xl, icode, bab1)

        bequ = dimo / (xl * xl * xl)
        if (icode == 1) then
            bdel = 1.0e-3_wp
            call me%findb0(stps, bdel, val, beq, rr0)
            if (val) bequ = beq
        end if
        bbx = norm2(b) / bequ

    end subroutine igrfc