get_sun_fraction Function

public function get_sun_fraction(b, rad_body, rad_sun, eph, et, rv, model, rbubble, use_geometric, info) result(phi)

Compute the "sun fraction" using the selected shadow model.

Arguments

Type IntentOptional Attributes Name
type(celestial_body), intent(in) :: b

eclipsing body

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

radius of the eclipsing body [km]

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

radius of the Sun [km]

class(ephemeris_class), intent(inout) :: eph

the ephemeris to use for sun and ssb (if necessary)

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

observer ephemeris time (sec)

real(kind=wp), intent(in), dimension(6) :: rv

state of the spacecraft (j2000-body frame)

integer, intent(in) :: model

algorithm to use:

  • 1: circular cubic shadow model
  • 2-4: solar fraction model
real(kind=wp), intent(in) :: rbubble

eclipse bubble [km]. see the reference. if rbubble=0, then no bubble is used. only used if model=1

logical, intent(in), optional :: use_geometric

if true, use geometric positions (no light time or stellar aberration correction) default = false

character(len=:), intent(out), optional, allocatable :: info

info string

Return Value real(kind=wp)

solar fraction returned:

  • if model=1, circular cubic sun frac value:
    • >0 no eclipse
    • <0 eclipse
    • =0 on the eclipse line
  • if model=2, true solar fraction value [0=total eclipse, 1=no eclipse], with model of umbra/penumbra/antumbra (Wertz, 1978)
  • if model=3, alternate version of solar fraction (Montenbruck and Gill)
  • if model=4, alternate version of solar fraction (nyxspace)

Calls

proc~~get_sun_fraction~~CallsGraph proc~get_sun_fraction get_sun_fraction get_r get_r proc~get_sun_fraction->get_r proc~apparent_position apparent_position proc~get_sun_fraction->proc~apparent_position proc~cubic_shadow_model cubic_shadow_model proc~get_sun_fraction->proc~cubic_shadow_model proc~from_j2000body_to_j2000ssb from_j2000body_to_j2000ssb proc~get_sun_fraction->proc~from_j2000body_to_j2000ssb proc~solar_fraction solar_fraction proc~get_sun_fraction->proc~solar_fraction proc~solar_fraction_alt solar_fraction_alt proc~get_sun_fraction->proc~solar_fraction_alt proc~solar_fraction_alt2 solar_fraction_alt2 proc~get_sun_fraction->proc~solar_fraction_alt2 proc~apparent_position->get_r proc~axis_angle_rotation axis_angle_rotation proc~apparent_position->proc~axis_angle_rotation proc~cross cross proc~apparent_position->proc~cross proc~unit unit proc~apparent_position->proc~unit acosd acosd proc~cubic_shadow_model->acosd asind asind proc~cubic_shadow_model->asind proc~cubic_shadow_model->proc~unit proc~wrap_angle wrap_angle proc~cubic_shadow_model->proc~wrap_angle proc~transform reference_frame%transform proc~from_j2000body_to_j2000ssb->proc~transform proc~solar_fraction->proc~unit proc~solar_fraction_alt->proc~unit proc~solar_fraction_alt2->proc~unit proc~axis_angle_rotation->proc~cross proc~axis_angle_rotation->proc~unit get_c_cdot get_c_cdot proc~transform->get_c_cdot proc~rvcto_rvcfrom_icrf rvcto_rvcfrom_icrf proc~transform->proc~rvcto_rvcfrom_icrf get_rv get_rv proc~rvcto_rvcfrom_icrf->get_rv proc~from_primary_to_center two_body_rotating_frame%from_primary_to_center proc~rvcto_rvcfrom_icrf->proc~from_primary_to_center proc~from_primary_to_center->get_rv

Source Code

    function get_sun_fraction(b, rad_body, rad_sun, eph, et, rv, model, rbubble, use_geometric, info) result (phi)

    type(celestial_body),intent(in)      :: b           !! eclipsing body
    real(wp),intent(in)                  :: rad_body    !! radius of the eclipsing body [km]
    real(wp),intent(in)                  :: rad_sun     !! radius of the Sun [km]
    class(ephemeris_class),intent(inout) :: eph         !! the ephemeris to use for sun and ssb (if necessary)
    real(wp),intent(in)                  :: et          !! observer ephemeris time (sec)
    real(wp),dimension(6),intent(in)     :: rv          !! state of the spacecraft (j2000-body frame)
    integer,intent(in)                   :: model       !! algorithm to use:
                                                        !!
                                                        !!  * 1: circular cubic shadow model
                                                        !!  * 2-4: solar fraction model
    real(wp),intent(in)                  :: rbubble     !! eclipse bubble [km]. see the reference.
                                                        !! if rbubble=0, then no bubble is used.
                                                        !! only used if model=1
    logical,intent(in),optional          :: use_geometric  !! if true, use geometric positions
                                                           !! (no light time or stellar aberration correction)
                                                           !! default = false
    character(len=:),allocatable,intent(out),optional :: info !! info string
    real(wp) :: phi !! solar fraction returned:
                    !!
                    !!  * if `model=1`, circular cubic sun frac value:
                    !!     * `>0` no eclipse
                    !!     * `<0` eclipse
                    !!     * `=0` on the eclipse line
                    !!  * if `model=2`, true solar fraction value [0=total eclipse, 1=no eclipse],
                    !!    with model of umbra/penumbra/antumbra (Wertz, 1978)
                    !!  * if `model=3`, alternate version of solar fraction (Montenbruck and Gill)
                    !!  * if `model=4`, alternate version of solar fraction (nyxspace)

    logical :: status_ok !! true if no problems
    real(wp),dimension(3) :: r_sun !! apparent state of the sun (j2000-ssb frame)
    real(wp),dimension(3) :: r_body !! apparent state of the eclipsing body (j2000-ssb frame)
    real(wp),dimension(6) :: rv_ssb !! state of the spacecraft !! (j2000-ssb frame)
    logical :: use_apparent

    if (present(use_geometric)) then
        use_apparent = .not. use_geometric
    else
        use_apparent = .true.
    end if

    if (use_apparent) then
        ! apparent position of sun and body wrt to the spacecraft
        call from_j2000body_to_j2000ssb(b, eph, et, rv, rv_ssb) ! state of spacecraft in j2000-ssb
        call apparent_position(eph, body_sun, et, rv_ssb, r_sun,  status_ok) ! apparent position of sun in j2000
        if (.not. status_ok) error stop 'error getting apparent sun position'
        call apparent_position(eph, b,        et, rv_ssb, r_body, status_ok) ! apparent position of body in j2000
        if (.not. status_ok) error stop 'error getting apparent body position'
    else
        ! use geometric positions
        r_body = -rv(1:3) ! geometric position of body wrt spacecraft in j2000
        call eph%get_r(et, body_sun, b, r_sun, status_ok) ! geometric position of sun wrt body in j2000
        if (.not. status_ok) error stop 'error getting geometric sun position'
        r_sun = r_body + r_sun ! geometric position of sun wrt spacecraft in j2000
    end if

    ! compute sun fraction value
    select case(model)
    case(1); call cubic_shadow_model( r_sun, rad_sun, r_body, rad_body, phi, rbubble)
    case(2); call solar_fraction(     r_sun, rad_sun, r_body, rad_body, phi, info)
    case(3); call solar_fraction_alt( r_sun, rad_sun, r_body, rad_body, phi, info)
    case(4); call solar_fraction_alt2(r_sun, rad_sun, r_body, rad_body, phi)
    case default
        error stop 'invalid sun fraction model'
    end select

    end function get_sun_fraction