Skip to content

Latest commit

 

History

History
175 lines (88 loc) · 5.1 KB

makegeoidgrid.pod

File metadata and controls

175 lines (88 loc) · 5.1 KB

MakeGeoidGrid

MakeGeoidGrid -

Create a global map of the geoid.

SYNOPSIS

SUBROUTINE MakeGeoidGrid (

GEOID, CILM, LMAX, R0POT, GM, POTREF, OMEGA, R, GRIDTYPE, ORDER, NLAT, NLONG, INTERVAL, LMAX_CALC, A, F )

    REAL*8

    GEOID(NLAT, NLONG), CILM(2, LMAX+1, LMAX+1), R0POT, GM, POTREF, OMEGA, R

    INTEGER

    LMAX, GRIDTYPE, ORDER, NLAT, NLONG

    REAL*8, OPTIONAL

    INTERVAL, A, F

    INTEGER, OPTIONAL

    LMAX_CALC

DESCRIPTION

MakeGeoidGrid will create a global map of the geoid, accurate to either first, second, or third order, using the method described in Wieczorek (2007; equation 19-20). The algorithm expands the potential in a Taylor series on a spherical interface of radius R, and computes the height above this interface to the potential POTREF exactly from the linear, quadratic, or cubic equation at each grid point. If the optional parameters A and F are specified, the geoid height will be referenced to a flattened ellipsoid with semi-major axis A and flattening F. The pseudo-rotational potential is explicitly accounted for by specifying the angular rotation rate OMEGA of the planet.

It should be noted that this geoid calculation is only strictly exact when the radius R lies above the maximum radius of the planet. Furthermore, the geoid is only strictly valid when it lies above the surface of the planet (it is necessary to know the density structure of the planet when calculating the potential below the surface).

The geoid can be computed on one of four different grids: (1) a Gauss-Legendre quadrature grid (see MakeGridGLQ), (2) A N by N equally sampled grid (see MakeGridDH), (3) an N by 2N equally spaced grid (see MakeGridDH), or (4) A 2D Cartesian grid (see MakeGrid2D),.

ARGUMENTS

GEOID

(output) REAL*8, DIMENSION(NLAT, NLONG)

A global grid of the height to the potential POTREF above a sphere of radius R (or above a flattened ellipsoid if both A and F are specified). The number of latitude and longitude points depends upon GRIDTYPE: (1) LMAX+1 by 2*LMAX + 1, (2) 2*LMAX+2 by 2*LMAX+2, (3) 2*LMAX+2 by 4*LMAX+4, or (4) 180/INTERVAL + 1 by 360/INTERVAL + 1.

CILM

(input) REAL*8, DIMENSION (2, LMAX+1, LMAX+1)

The real spherical harmonic coefficients (geodesy normalized) of the gravitational potential referenced to a spherical interface of radius R0POT.

LMAX

(input) INTEGER

The maximum spherical harmonic degree of the gravitational-potential coefficients. For GRIDTYPEs 1, 2 and 3, this determines the number of latitudinal and longitudinal samples.

R0POT

(input) REAL*8

The reference radius of the spherical harmonic coefficients.

GM

(input) REAL*8

The product of the gravitational constant and mass of the planet.

POTREF

(input) REAL*8

The value of the potential on the chosen geoid, in SI units.

OMEGA

(input) REAL*8

The angular rotation rate of the planet.

R

(input) REAL*8

The radius of the reference sphere that the Taylor expansion of the potential is performed on. If A and F are not specified, the geoid height will be referenced to this spherical interface.

GRIDTYPE

(input) INTEGER

The output grid is (1) a Gauss-Legendre quadrature grid whose grid nodes are determined by LMAX, (2) an equally sampled N by N grid used with the Driscoll and Healy (1994) sampling theorem, (3) ar a similar N by 2N grid that is oversampled in longitude, or (4) a 2D Cartesian grid with latitudinal and longitudinal spacing given by INTERVAL.

ORDER

(input) INTEGER

The order of the Taylor series expansion of the potential about the reference radius R. This can be either 1, 2, or 3.

NLAT

(output) INTEGER

The number of latitudinal samples.

NLONG

(output) INTEGER

The number of longitudinal samples.

INTERVAL

(input) REAL*8, OPTIONAL

The latitudinal and longitudinal spacing of the output grid for GRIDTYPE is 4.

LMAX_CALC

(input) INTEGER, OPTIONAL

The maximum degree used in evaluating the spherical harmonic coefficients. This must be less than or equal to LMAX.

A

(input) REAL*8, OPTIONAL

The semi-major axis of the flattened ellipsoid that the output grid GEOID is referenced to. The optional parameter F must also be specified.

F

(input) REAL*8, OPTIONAL

The flattening (R_equator - R_pole)/R_equator of the reference ellipsoid. The optional parameter A (i.e., R_equator) must be specified.

NOTES

This routine uses geodesy 4-pi normalized spherical harmonics that exclude the Condon-Shortley phase. This can not be modified.

SEE ALSO

makegrid2d(1), makegridglq(1), makegriddh(1)

http://shtools.ipgp.fr/

REFERENCES

Driscoll, J.R. and D.M. Healy, Computing Fourier transforms and convolutions on the 2-sphere, Adv. Appl. Math., 15, 202-250, 1994.

Wieczorek, M. A. Gravity and topography of the terrestrial planets, Treatise on Geophysics, 10, 165-206, 2007.

COPYRIGHT AND LICENSE

Copyright 2012 by Mark Wieczorek <wieczor@ipgp.fr>.

This is free software; you can distribute and modify it under the terms of the revised BSD license.