getshc Subroutine

private subroutine getshc(Fspec, Nmax, Erad, Gh, Ier)

Reads spherical harmonic coefficients from the specified file into an array.

Author

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

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: Fspec

File specification

integer, intent(out) :: Nmax

Maximum degree and order of model

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

Earth's radius associated with the spherical harmonic coefficients, in the same units as elevation

real(kind=wp), intent(out), dimension(*) :: Gh

Schmidt quasi-normal internal spherical harmonic coefficients

integer, intent(out) :: Ier

Error number:

  • 0, no error
  • -2, records out of order
  • FORTRAN run-time error number

Called by

proc~~getshc~~CalledByGraph proc~getshc shellig_module::getshc proc~feldcof shellig_module::shellig_type%feldcof proc~feldcof->proc~getshc 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 getshc(Fspec, Nmax, Erad, Gh, Ier)

        character(len=*), intent(in) :: Fspec !! File specification
        integer, intent(out) :: Nmax !! Maximum degree and order of model
        real(wp), intent(out) :: Erad !! Earth's radius associated with the spherical
                                !! harmonic coefficients, in the same units as
                                !! elevation
        real(wp), dimension(*), intent(out) :: Gh !! Schmidt quasi-normal internal spherical
                                           !! harmonic coefficients
        integer, intent(out) :: Ier !! Error number:
                              !!
                              !!  * 0, no error
                              !!  * -2, records out of order
                              !!  * FORTRAN run-time error number

        integer :: iu !! logical unit number
        real(wp) :: g, h
        integer :: i, m, mm, n, nn

        read_file: block
            ! ---------------------------------------------------------------
            !  Open coefficient file. Read past first header record.
            !  Read degree and order of model and Earth's radius.
            ! ---------------------------------------------------------------
            open (newunit=Iu, FILE=Fspec, STATUS='OLD', IOSTAT=Ier)
            if (Ier /= 0) then
                write (*, *) 'Error opening file: '//trim(fspec)
                exit read_file
            end if
            read (Iu, *, IOSTAT=Ier)
            if (Ier /= 0) exit read_file
            read (Iu, *, IOSTAT=Ier) Nmax, Erad
            if (Ier /= 0) exit read_file

            ! ---------------------------------------------------------------
            !  Read the coefficient file, arranged as follows:
            !
            !          N     M     G     H
            !          ----------------------
            !            /   1     0    GH(1)  -
            !           /  1     1    GH(2) GH(3)
            !          /  2     0    GH(4)  -
            !         /  2     1    GH(5) GH(6)
            !      NMAX*(NMAX+3)/2   /  2     2    GH(7) GH(8)
            !         records    \  3     0    GH(9)  -
            !         \      .     .     .     .
            !          \  .     .     .     .
            !      NMAX*(NMAX+2)     \  .     .     .     .
            !      elements in GH      \  NMAX  NMAX   .     .
            !
            !  N and M are, respectively, the degree and order of the
            !  coefficient.
            ! ---------------------------------------------------------------
            i = 0
            main: do nn = 1, Nmax
                do mm = 0, nn
                    read (Iu, *, IOSTAT=Ier) n, m, g, h
                    if (Ier /= 0) exit main
                    if (nn /= n .or. mm /= m) then
                        Ier = -2
                        exit main
                    end if
                    i = i + 1
                    Gh(i) = g
                    if (m /= 0) then
                        i = i + 1
                        Gh(i) = h
                    end if
                end do
            end do main

        end block read_file

        close (Iu)

    end subroutine getshc