shellig_module Module

IGRF model

History

  • SHELLIG.FOR, Version 2.0, January 1992
  • 11/01/91-DKB- SHELLG: lowest starting point for B0 search is 2
  • 1/27/92-DKB- Adopted to IGRF-91 coefficients model
  • 2/05/92-DKB- Reduce variable-names: INTER(P)SHC,EXTRA(P)SHC,INITI(ALI)ZE
  • 8/08/95-DKB- Updated to IGRF-45-95; new coeff. DGRF90, IGRF95, IGRF95S
  • 5/31/00-DKB- Updated to IGRF-45-00; new coeff.: IGRF00, IGRF00s
  • 3/24/05-DKB- Updated to IGRF-45-10; new coeff.: IGRF05, IGRF05s

Uses

  • module~~shellig_module~~UsesGraph module~shellig_module shellig_module module~radbelt_kinds_module radbelt_kinds_module module~shellig_module->module~radbelt_kinds_module iso_fortran_env iso_fortran_env module~radbelt_kinds_module->iso_fortran_env

Used by

  • module~~shellig_module~~UsedByGraph module~shellig_module shellig_module module~radbelt_module radbelt_module module~radbelt_module->module~shellig_module module~radbelt_c_module radbelt_c_module module~radbelt_c_module->module~radbelt_module

Variables

Type Visibility Attributes Name Initial
integer, private, parameter :: filename_len = 14

length of the model data file names

real(kind=wp), private, parameter :: Era = 6371.2_wp

earth radius for normalization of cartesian coordinates (6371.2 km)

real(kind=wp), private, parameter :: erequ = 6378.16_wp
real(kind=wp), private, parameter :: erpol = 6356.775_wp
real(kind=wp), private, parameter :: Aquad = erequ*erequ

square of major half axis for earth ellipsoid as recommended by international astronomical union

real(kind=wp), private, parameter :: Bquad = erpol*erpol

square of minor half axis for earth ellipsoid as recommended by international astronomical union

real(kind=wp), private, parameter :: Umr = atan(1.0_wp)*4.0_wp/180.0_wp

atan(1.0)4./180. umr=

real(kind=wp), private, parameter, dimension(3, 3) :: u = reshape([+0.3511737_wp, -0.9148385_wp, -0.1993679_wp, +0.9335804_wp, +0.3583680_wp, +0.0000000_wp, +0.0714471_wp, -0.1861260_wp, +0.9799247_wp], [3, 3])
integer, private, parameter :: max_loop_index = 3333

used in shellg for the field line integration loop


Derived Types

type, public ::  shellig_type

Components

Type Visibility Attributes Name Initial
character(len=:), private, allocatable :: igrf_dir

directory containing the data files

real(kind=wp), private, dimension(3) :: sp = 0.0_wp
real(kind=wp), private, dimension(3) :: xi = 0.0_wp
real(kind=wp), private, dimension(144) :: h = 0.0_wp

Field model coefficients adjusted for shellg

integer, private :: iyea = 0

the int year corresponding to the file name that has been read

character(len=:), private, allocatable :: name

file name

integer, private :: nmax = 0

maximum order of spherical harmonics

real(kind=wp), private :: Time = 0.0_wp

year (decimal: 1973.5) for which magnetic field is to be calculated

real(kind=wp), private, dimension(144) :: g = 0.0_wp

g(m) -- normalized field coefficients (see feldcof) m=nmax*(nmax+2)

integer, private :: nmax1 = 0

saved variables from the file

integer, private :: nmax2 = 0

saved variables from the file

real(kind=wp), private, dimension(144) :: g_cache = 0.0_wp

saved g from the file

real(kind=wp), private :: step = 0.20_wp

step size for field line tracing

real(kind=wp), private :: steq = 0.03_wp

step size for integration

real(kind=wp), private, dimension(120) :: gh2 = 0.0_wp
real(kind=wp), private, dimension(:, :), allocatable :: p

this was p(8,100) in the original code. used for the field line integration loop. changed it to be allocatable since it was changed to be p(8,3334).

Type-Bound Procedures

procedure, public :: igrfc
procedure, public :: igrf
procedure, public :: feldcof
procedure, public :: feldc
procedure, public :: feldg
procedure, public :: shellc
procedure, public :: shellg
procedure, public :: findb0
procedure, private :: feldi
procedure, private :: stoer
procedure, public :: get_data_file_dir
procedure, public :: set_data_file_dir
procedure, public :: destroy => destroy_shellig_type

Functions

private function get_data_file_dir(me) result(dir)

Get the directory containing the data files.

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(in) :: me

Return Value character(len=:), allocatable

private pure function geo_to_cart(glat, glon, alt) result(x)

geodetic to scaled cartesian coordinates

Arguments

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

geodetic latitude in degrees (north)

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

geodetic longitude in degrees (east)

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

altitude in km above sea level

Return Value real(kind=wp), dimension(3)

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

Subroutines

private subroutine destroy_shellig_type(me)

Destroy a shellig_type.

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(out) :: me

private subroutine set_data_file_dir(me, dir)

Set the directory containing the data files.

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(inout) :: me
character(len=*), intent(in) :: dir

private subroutine igrf(me, lon, lat, height, year, xl, bbx)

Wrapper for IGRF functions.

Arguments

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

geodetic longitude in degrees (east)

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

geodetic latitude in degrees (north)

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

altitude in km above sea level

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

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

Alternate version of igrf for cartesian coordinates.

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

private subroutine findb0(me, stps, bdel, value, bequ, rr0)

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(inout) :: me
real(kind=wp), intent(in) :: stps
real(kind=wp), intent(inout) :: bdel
logical, intent(out) :: value
real(kind=wp), intent(out) :: bequ
real(kind=wp), intent(out) :: rr0

private subroutine shellc(me, v, dimo, fl, icode, b0)

Wrapper to shellg to be used with cartesian coordinates.

Read more…

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) :: dimo

dipol moment in gauss (normalized to earth radius)

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

l-value

integer, intent(out) :: icode Read more…
real(kind=wp), intent(out) :: b0

magnetic field strength in gauss

private subroutine shellg(me, glat, glon, alt, dimo, fl, icode, b0, v)

calculates l-value for specified geodaetic coordinates, altitude and gemagnetic field model.

Read more…

Arguments

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

geodetic latitude in degrees (north)

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

geodetic longitude in degrees (east)

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

altitude in km above sea level

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

dipol moment in gauss (normalized to earth radius)

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

l-value

integer, intent(out) :: icode Read more…
real(kind=wp), intent(out) :: b0

magnetic field strength in gauss

real(kind=wp), intent(in), optional, dimension(3) :: v

cartesian coordinates in earth radii (6371.2 km)

Read more…

private subroutine stoer(me, p, bq, r)

subroutine used for field line tracing in shellg. calls entry point feldi in geomagnetic field subroutine feldg

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(inout) :: me
real(kind=wp), intent(inout), dimension(7) :: p
real(kind=wp), intent(out) :: bq
real(kind=wp), intent(out) :: r

private subroutine feldg(me, glat, glon, alt, bnorth, beast, bdown, Babs)

Calculates earth magnetic field from spherical harmonics model

Read more…

Arguments

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

geodetic latitude in degrees (north)

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

geodetic longitude in degrees (east)

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

altitude in km above sea level

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

components of the field with respect to the local geodetic coordinate system, with axis pointing in the tangential plane to the north, east and downward.

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

components of the field with respect to the local geodetic coordinate system, with axis pointing in the tangential plane to the north, east and downward.

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

components of the field with respect to the local geodetic coordinate system, with axis pointing in the tangential plane to the north, east and downward.

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

magnetic field strength in gauss

private subroutine feldc(me, v, b)

Alternate version of feldg to be used with cartesian coordinates

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

field components

private subroutine feldi(me)

Used for l computation.

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(inout) :: me

private subroutine feldcof(me, year, dimo)

Determines coefficients and dipol moment from IGRF models

Read more…

Arguments

Type IntentOptional Attributes Name
class(shellig_type), intent(inout) :: me
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) :: dimo

geomagnetic dipol moment in gauss (normalized to earth's radius) at the time (year)

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

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

Read more…

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:

Read more…

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

Interpolates linearly, in time, between two spherical harmonic models.

Read more…

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 earlier model

integer, intent(in) :: nmax1

Maximum degree and order of earlier model

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

Schmidt quasi-normal internal spherical harmonic coefficients of earlier model

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

Date of later model

integer, intent(in) :: nmax2

Maximum degree and order of later model

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

Schmidt quasi-normal internal spherical harmonic coefficients of later model

integer, intent(out) :: nmax

Maximum degree and order of resulting model

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

Coefficients of resulting model

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

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

Read more…

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