extrashc Subroutine

private subroutine extrashc(date, dte1, nmax1, gh1, nmax2, gh2, nmax, gh)

Extrapolates linearly a spherical harmonic model with a rate-of-change model.

The coefficients (GH) of the resulting model, at date DATE, are computed by linearly extrapolating the coef- ficients of the base model (GH1), at date DTE1, using those of the rate-of-change model (GH2), at date DTE2. If one model is smaller than the other, the extrapolation is performed with the missing coefficients assumed to be 0.

Author

  • Version 1.01, A. Zunde USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: date

Date of resulting model (in decimal year)

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

Date of base model

integer, intent(in) :: nmax1

Maximum degree and order of base model

real(kind=wp), intent(in) :: gh1(*)

Schmidt quasi-normal internal spherical harmonic coefficients of base model

integer, intent(in) :: nmax2

Maximum degree and order of rate-of-change model

real(kind=wp), intent(in) :: gh2(*)

Schmidt quasi-normal internal spherical harmonic coefficients of rate-of-change model

integer, intent(out) :: nmax

Maximum degree and order of resulting model

real(kind=wp), intent(out) :: gh(*)

Coefficients of resulting model


Called by

proc~~extrashc~~CalledByGraph proc~extrashc shellig_module::extrashc proc~feldcof shellig_module::shellig_type%feldcof proc~feldcof->proc~extrashc proc~igrf shellig_module::shellig_type%igrf proc~igrf->proc~feldcof proc~igrfc shellig_module::shellig_type%igrfc proc~igrfc->proc~feldcof proc~get_flux_c_ radbelt_module::radbelt_type%get_flux_c_ proc~get_flux_c_->proc~igrfc proc~get_flux_g_ radbelt_module::radbelt_type%get_flux_g_ proc~get_flux_g_->proc~igrf none~get_flux radbelt_module::radbelt_type%get_flux none~get_flux->proc~get_flux_c_ none~get_flux->proc~get_flux_g_ 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 extrashc(date, dte1, nmax1, gh1, nmax2, gh2, nmax, gh)

        real(wp), intent(in) :: date !! Date of resulting model (in decimal year)
        real(wp), intent(in) :: dte1 !! Date of base model
        integer, intent(in) :: nmax1 !! Maximum degree and order of base model
        real(wp), intent(in) :: gh1(*) !! Schmidt quasi-normal internal spherical harmonic coefficients of base model
        integer, intent(in) :: nmax2 !! Maximum degree and order of rate-of-change model
        real(wp), intent(in) :: gh2(*) !! Schmidt quasi-normal internal spherical harmonic coefficients of rate-of-change model
        real(wp), intent(out) :: gh(*) !! Coefficients of resulting model
        integer, intent(out) :: nmax !! Maximum degree and order of resulting model

        real(wp) :: factor
        integer :: i, k, l

        factor = (date - dte1)

        if (nmax1 == nmax2) then
            k = nmax1 * (nmax1 + 2)
            nmax = nmax1
        elseif (nmax1 > nmax2) then
            k = nmax2 * (nmax2 + 2)
            l = nmax1 * (nmax1 + 2)
            do i = k + 1, l
                gh(i) = gh1(i)
            end do
            nmax = nmax1
        else
            k = nmax1 * (nmax1 + 2)
            l = nmax2 * (nmax2 + 2)
            do i = k + 1, l
                gh(i) = factor * gh2(i)
            end do
            nmax = nmax2
        end if

        do i = 1, k
            gh(i) = gh1(i) + factor * gh2(i)
        end do

    end subroutine extrashc