grid_tools  1.7.0
 All Data Structures Files Functions Variables
pmat5 Module Reference

Utility routines for orienting the globe and basic geographical mappings. More...

Data Types

interface  ctoc_schm
 
interface  ctog
 
interface  ctogr
 
interface  frametwist
 
interface  grtoc
 
interface  gtoc
 
interface  gtoframe
 
interface  ininmap
 
interface  inivmap
 
interface  paraframe
 
interface  plctoc
 
interface  plrot
 
interface  plroti
 

Public Member Functions

subroutine dctoc (s, xc1, xc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s. More...
 
subroutine dctocd (s, xc1, xc2, dxc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, and its jacobian, dxc2. More...
 
subroutine dctocdd (s, xc1, xc2, dxc2, ddxc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, its jacobian, dxc2, and its 2nd derivative, ddxc2. More...
 
subroutine dctog (xe, dlat, dlon)
 Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine dctogr (xe, rlat, rlon)
 Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (radians) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine dgtoc (dlat, dlon, xe)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine dininmap (alon0, alat0, rot3)
 Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) on the projection axis. More...
 
subroutine dinivmap (alon0, alat0, rot3)
 Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) at the viewing nadir. More...
 
subroutine, public frametwist (xp, yp, zp, xv, yv, zv, twist)
 Given a principal cartesian orthonormal frame, {xp,yp,zp} (i.e., at P with Earth-centered cartesians, zp), and another similar frame {xv,yv,zv} at V with Earth-centered cartesians zv, find the relative rotation angle, "twist" by which the frame at V is rotated in the counterclockwise sense relative to the parallel-transportation of P's frame to V. More...
 
subroutine, public paraframe (xp, yp, zp, xv, yv, zv)
 Take a principal reference orthonormal frame, {xp,yp,zp} and a dependent point defined by unit vector, zv, and complete the V-frame cartesian components, {xv,yv}, that are the result of parallel-transport of {xp,yp} along the geodesic between P and V. More...
 
subroutine, public plctoc (s, n, x, y, z)
 Perform Schmidt transformation with scaling parameter s to a polyline. More...
 
subroutine, public plrot (rot3, n, x, y, z)
 Apply a constant rotation to a three dimensional polyline. More...
 
subroutine, public plroti (rot3, n, x, y, z)
 Invert the rotation of a three-dimensional polyline. More...
 
subroutine sctocd (s, xc1, xc2, dxc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, and its jacobian, dxc2. More...
 
subroutine sctocdd (s, xc1, xc2, dxc2, ddxc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, its jacobian, dxc2, and its 2nd derivative, ddxc2. More...
 
subroutine sctog (xe, dlat, dlon)
 Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine sctogr (xe, rlat, rlon)
 Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (radians) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine sframetwist (sxp, syp, szp, sxv, syv, szv, stwist)
 Given a principal cartesian orthonormal frame, {xp,yp,zp} (i.e., at P with Earth-centered cartesians, zp), and another similar frame {xv,yv,zv} at V with Earth-centered cartesians zv, find the relative rotation angle, "twist" by which the frame at V is rotated in the counterclockwise sense relative to the parallel-transportation of P's frame to V. More...
 
subroutine sgrtoc (rlat, rlon, xe)
 Transform "Geographical" to "Cartesian" coordinates. More...
 
subroutine sgrtocdd (rlat, rlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
 Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivatives of cartesians wrt latitude and longitude. More...
 
subroutine sgtoc (dlat, dlon, xe)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine sgtocd (dlat, dlon, xe, dxedlat, dxedlon)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine sgtocdd (dlat, dlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine sgtoframem (splat, splon, sorth)
 From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location as an orthonormal matrix, orth. More...
 
subroutine sinivmap (alon0, alat0, rot3)
 Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) at the viewing nadir. More...
 
subroutine sparaframe (sxp, syp, szp, sxv, syv, szv)
 Take a principal reference orthonormal frame, {xp,yp,zp} and a dependent point defined by unit vector, zv, and complete the V-frame cartesian components, {xv,yv}, that are the result of parallel-transport of {xp,yp} along the geodesic between P and V. More...
 

Private Member Functions

subroutine dgrtoc (rlat, rlon, xe)
 Transform "Geographical" to "Cartesian" coordinates. More...
 
subroutine dgrtocd (rlat, rlon, xe, dxedlat, dxedlon)
 Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude. More...
 
subroutine dgrtocdd (rlat, rlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
 Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivatives of cartesians wrt latitude and longitude. More...
 
subroutine dgtocd (dlat, dlon, xe, dxedlat, dxedlon)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine dgtocdd (dlat, dlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
 Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian. More...
 
subroutine dplrot (rot3, n, x, y, z)
 Apply a constant rotation to a three dimensional polyline. More...
 
subroutine dplroti (rot3, n, x, y, z)
 Invert the rotation of a three-dimensional polyline. More...
 
subroutine gtoframem (plat, plon, orth)
 From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location as an orthonormal matrix, orth. More...
 
subroutine gtoframev (plat, plon, xp, yp, zp)
 Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal cartesian frame, {xp,yp,zp} in earth-centered coordinates. More...
 
subroutine sctoc (s, xc1, xc2)
 Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s. More...
 
subroutine sgrtocd (rlat, rlon, xe, dxedlat, dxedlon)
 Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude. More...
 
subroutine sgtoframev (splat, splon, sxp, syp, szp)
 Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal cartesian frame, {xp,yp,zp} in earth-centered coordinates. More...
 
subroutine sininmap (alon0, alat0, rot3)
 Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) on the projection axis. More...
 

Detailed Description

Utility routines for orienting the globe and basic geographical mappings.

Author
R. J. Purser

Definition at line 27 of file pmat5.f90.

Member Function/Subroutine Documentation

subroutine pmat5::dctoc ( real(dp), intent(in)  s,
real(dp), dimension(3), intent(inout)  xc1,
real(dp), dimension(3), intent(inout)  xc2 
)

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s.

Parameters
[in]sSchmidt scaling parameter.
[in,out]xc1input cartesian 3-vector.
[in,out]xc2output cartesian unit 3-vector.
Author
R. J. Purser

Definition at line 915 of file pmat5.f90.

subroutine pmat5::dctocd ( real(dp), intent(in)  s,
real(dp), dimension(3), intent(inout)  xc1,
real(dp), dimension(3), intent(inout)  xc2,
real(dp), dimension(3,3), intent(out)  dxc2 
)

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, and its jacobian, dxc2.

Parameters
[in]sSchmidt scaling parameter.
[in,out]xc1input cartesian 3-vector.
[in,out]xc2output cartesian unit 3-vector.
[out]dxc2jacobian 1st derivative, d(xc2)/d(xc1).
Author
R. J. Purser

Definition at line 946 of file pmat5.f90.

subroutine pmat5::dctocdd ( real(dp), intent(in)  s,
real(dp), dimension(3), intent(inout)  xc1,
real(dp), dimension(3), intent(inout)  xc2,
real(dp), dimension(3,3), intent(out)  dxc2,
real(dp), dimension(3,3,3), intent(out)  ddxc2 
)

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, its jacobian, dxc2, and its 2nd derivative, ddxc2.

Parameters
[in]sSchmidt scaling parameter.
[in,out]xc1input cartesian 3-vector.
[in,out]xc2output cartesian unit 3-vector.
[out]dxc2jacobian 1st derivative, d(xc2)/d(xc1).
[out]ddxc22nd derivative, d^2(xc2)/(d(xc1)d(xc1)).
Author
R. J. Purser

Definition at line 987 of file pmat5.f90.

subroutine pmat5::dctog ( real(dp), dimension(3), intent(in)  xe,
real(dp), intent(out)  dlat,
real(dp), intent(out)  dlon 
)

Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Double precision version.

Parameters
[in]xeEarth-centered cartesian unit 3-vector.
[out]dlatdegrees latitude.
[out]dlondegrees longitude.
Author
R. J. Purser

Definition at line 397 of file pmat5.f90.

subroutine pmat5::dctogr ( real(dp), dimension(3), intent(in)  xe,
real(dp), intent(out)  rlat,
real(dp), intent(out)  rlon 
)

Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (radians) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Double precision version.

Parameters
[in]xeEarth-centered cartesian unit 3-vector.
[out]rlatradians latitude.
[out]rlonradians longitude.
Author
R. J. Purser

Definition at line 197 of file pmat5.f90.

subroutine pmat5::dgrtoc ( real(dp), intent(in)  rlat,
real(dp), intent(in)  rlon,
real(dp), dimension(3), intent(out)  xe 
)
private

Transform "Geographical" to "Cartesian" coordinates.

Double proecision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
Author
R. J. Purser

Definition at line 236 of file pmat5.f90.

subroutine pmat5::dgrtocd ( real(dp), intent(in)  rlat,
real(dp), intent(in)  rlon,
real(dp), dimension(3), intent(out)  xe,
real(dp), dimension(3), intent(out)  dxedlat,
real(dp), dimension(3), intent(out)  dxedlon 
)
private

Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude.

Double precision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
[out]dxedlatDerivative, d(xe)/d(rlat).
[out]dxedlonDerivative, d(xe)/d(rlon).
Author
R. J. Purser

Definition at line 279 of file pmat5.f90.

subroutine pmat5::dgrtocdd ( real(dp), intent(in)  rlat,
real(dp), intent(in)  rlon,
real(dp), dimension(3), intent(out)  xe,
real(dp), dimension(3), intent(out)  dxedlat,
real(dp), dimension(3), intent(out)  dxedlon,
real(dp), dimension(3), intent(out)  ddxedlatdlat,
real(dp), dimension(3), intent(out)  ddxedlatdlon,
real(dp), dimension(3), intent(out)  ddxedlondlon 
)
private

Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivatives of cartesians wrt latitude and longitude.

Double precision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
[out]dxedlatDerivative, d(xe)/d(rlat).
[out]dxedlonDerivative, d(xe)/d(rlon).
[out]ddxedlatdlatDerivative, d^2(xe)/(d(rlat)d(rlat)).
[out]ddxedlatdlonDerivative, d^2(xe)/(d(rlat)d(rlon)).
[out]ddxedlondlonDerivative, d^2(xe)/(d(rlon)d(rlon)).
Author
R. J. Purser

Definition at line 338 of file pmat5.f90.

subroutine pmat5::dgtoc ( real(dp), intent(in)  dlat,
real(dp), intent(in)  dlon,
real(dp), dimension(3), intent(out)  xe 
)

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Double precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xecartesian unit 3-vector.
Author
R. J. Purser

Definition at line 444 of file pmat5.f90.

subroutine pmat5::dgtocd ( real(dp), intent(in)  dlat,
real(dp), intent(in)  dlon,
real(dp), dimension(3), intent(out)  xe,
real(dp), dimension(3), intent(out)  dxedlat,
real(dp), dimension(3), intent(out)  dxedlon 
)
private

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Also, return the partial derivatives of xe wrt latitude and longitude. Double precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xecartesian unit 3-vector.
[out]dxedlatderivative, d(xe)/d(dlat).
[out]dxedlonderivative, d(xe)/d(dlon).
Author
R. J. Purser

Definition at line 495 of file pmat5.f90.

subroutine pmat5::dgtocdd ( real(dp), intent(in)  dlat,
real(dp), intent(in)  dlon,
real(dp), dimension(3), intent(out)  xe,
real(dp), dimension(3), intent(out)  dxedlat,
real(dp), dimension(3), intent(out)  dxedlon,
real(dp), dimension(3), intent(out)  ddxedlatdlat,
real(dp), dimension(3), intent(out)  ddxedlatdlon,
real(dp), dimension(3), intent(out)  ddxedlondlon 
)
private

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Also, return the 1st and 2nd partial derivatives of xe wrt latitude and longitude. Double precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xecartesian unit 3-vector.
[out]dxedlatd(xe)/d(dlat).
[out]dxedlond(xe)/d(dlon).
[out]ddxedlatdlatderivative, d^2(xe)/(d(dlat)d(dlat)).
[out]ddxedlatdlonderivative, d^2(xe)/(d(dlat)d(dlon)).
[out]ddxedlondlonderivative, d^2(xe)/(d(dlon)d(dlon)).
Author
R. J. Purser

Definition at line 563 of file pmat5.f90.

subroutine pmat5::dininmap ( real(dp), intent(in)  alon0,
real(dp), intent(in)  alat0,
real(dp), dimension(3,3), intent(out)  rot3 
)

Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) on the projection axis.

Double precision version.

Parameters
[in]alon0geographical longitude (degrees) of projection center.
[in]alat0geographical latitude (degrees) of projection center.
[out]rot3rotation matrix.
Author
R. J. Purser

Definition at line 85 of file pmat5.f90.

subroutine pmat5::dinivmap ( real(dp), intent(in)  alon0,
real(dp), intent(in)  alat0,
real(dp), dimension(3,3), intent(out)  rot3 
)

Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) at the viewing nadir.

Double precision version.

Parameters
[in]alon0geographical longitude (degrees) of viewing nadir.
[in]alat0geographical latitude (degrees) of viewing nadir.
[out]rot3rotation matrix.
Author
R. J. Purser

Definition at line 139 of file pmat5.f90.

subroutine pmat5::dplrot ( real(dp), dimension(3,3), intent(in)  rot3,
integer, intent(in)  n,
real(dp), dimension(n), intent(inout)  x,
real(dp), dimension(n), intent(inout)  y,
real(dp), dimension(n), intent(inout)  z 
)
private

Apply a constant rotation to a three dimensional polyline.

Parameters
[in]rot3rotation matrix.
[in]nnumber of points in the polyline.
[in,out]xcartesian components of the three dimensional polyline.
[in,out]ycartesian components of the three dimensional polyline.
[in,out]zcartesian components of the three dimensional polyline.
Author
R. J. Purser

Definition at line 1083 of file pmat5.f90.

subroutine pmat5::dplroti ( real(dp), dimension(3,3), intent(in)  rot3,
integer, intent(in)  n,
real(dp), dimension(n), intent(inout)  x,
real(dp), dimension(n), intent(inout)  y,
real(dp), dimension(n), intent(inout)  z 
)
private

Invert the rotation of a three-dimensional polyline.

Parameters
[in]rot3rotation to be inverted.
[in]nnumber of points in the polyline.
[in,out]xcartesian components of the three dimensional polyline.
[in,out]ycartesian components of the three dimensional polyline.
[in,out]zcartesian components of the three dimensional polyline.
Author
R. J. Purser

Definition at line 1105 of file pmat5.f90.

subroutine, public pmat5::frametwist ( real(dp), dimension(3), intent(in)  xp,
real(dp), dimension(3), intent(in)  yp,
real(dp), dimension(3), intent(in)  zp,
real(dp), dimension(3), intent(in)  xv,
real(dp), dimension(3), intent(in)  yv,
real(dp), dimension(3), intent(in)  zv,
real(dp), intent(out)  twist 
)

Given a principal cartesian orthonormal frame, {xp,yp,zp} (i.e., at P with Earth-centered cartesians, zp), and another similar frame {xv,yv,zv} at V with Earth-centered cartesians zv, find the relative rotation angle, "twist" by which the frame at V is rotated in the counterclockwise sense relative to the parallel-transportation of P's frame to V.

Note that, by symmetry, transposing P and V leads to the opposite twist. Double precision version.

Parameters
[in]xpP-frame cartesian X-vector.
[in]ypP-frame cartesian Y-vector.
[in]zpP-frame cartesian Z-vector.
[in]xvV-frame cartesian X-vector.
[in]yvV-frame cartesian Y-vector.
[in]zvV-frame cartesian Z-vector.
[out]twistrelative rotation angle (radians) of frames.
Author
R. J. Purser

Definition at line 777 of file pmat5.f90.

subroutine pmat5::gtoframem ( real(dp), intent(in)  plat,
real(dp), intent(in)  plon,
real(dp), dimension(3,3), intent(out)  orth 
)
private

From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location as an orthonormal matrix, orth.

Double precision version.

Parameters
[in]platlatitude (degrees) of point.
[in]plonlongitude (degrees) of point.
[out]orthorthonormal matrix.
Author
R. J. Purser

Definition at line 616 of file pmat5.f90.

References pmat5::gtoframe::gtoframev().

subroutine pmat5::gtoframev ( real(dp), intent(in)  plat,
real(dp), intent(in)  plon,
real(dp), dimension(3), intent(out)  xp,
real(dp), dimension(3), intent(out)  yp,
real(dp), dimension(3), intent(out)  zp 
)
private

Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal cartesian frame, {xp,yp,zp} in earth-centered coordinates.

Double precision version.

Parameters
[in]platlatitude (degrees) of point.
[in]plonlongitude (degrees) of point.
[out]xpunit X-basis vector of cartesian frame.
[out]ypunit Y-basis vector of cartesian frame.
[out]zpunit Z-basis vector of cartesian frame.
Author
R. J. Purser

Definition at line 661 of file pmat5.f90.

subroutine, public pmat5::paraframe ( real(dp), dimension(3), intent(in)  xp,
real(dp), dimension(3), intent(in)  yp,
real(dp), dimension(3), intent(in)  zp,
real(dp), dimension(3), intent(out)  xv,
real(dp), dimension(3), intent(out)  yv,
real(dp), dimension(3), intent(in)  zv 
)

Take a principal reference orthonormal frame, {xp,yp,zp} and a dependent point defined by unit vector, zv, and complete the V-frame cartesian components, {xv,yv}, that are the result of parallel-transport of {xp,yp} along the geodesic between P and V.

Double precision version.

Parameters
[in]xpreference orthonormal P-frame cartesian X-vector.
[in]ypreference orthonormal P-frame cartesian Y-vector.
[in]zpreference orthonormal P-frame cartesian Z-vector.
[out]xvV-frame cartesian X-vector.
[out]yvV-frame cartesian Y-vector.
[in]zvdependent point zenith, V-frame cartesian Z-vector.
Author
R. J. Purser

Definition at line 718 of file pmat5.f90.

subroutine, public pmat5::plctoc ( real(sp), intent(in)  s,
integer, intent(in)  n,
real(sp), dimension(n), intent(inout)  x,
real(sp), dimension(n), intent(inout)  y,
real(sp), dimension(n), intent(inout)  z 
)

Perform Schmidt transformation with scaling parameter s to a polyline.

Parameters
[in]sSchmidt scaling parameter.
[in]nnumber of points in the polyline.
[in,out]xcartesian components of the three dimensional polyline.
[in,out]ycartesian components of the three dimensional polyline.
[in,out]zcartesian components of the three dimensional polyline.
Author
R. J. Purser

Definition at line 1127 of file pmat5.f90.

subroutine, public pmat5::plrot ( real(sp), dimension(3,3), intent(in)  rot3,
integer, intent(in)  n,
real(sp), dimension(n), intent(inout)  x,
real(sp), dimension(n), intent(inout)  y,
real(sp), dimension(n), intent(inout)  z 
)

Apply a constant rotation to a three dimensional polyline.

Parameters
[in]rot3rotation matrix.
[in]nnumber of points in the polyline.
[in,out]xcartesian components of the three dimensional polyline.
[in,out]ycartesian components of the three dimensional polyline.
[in,out]zcartesian components of the three dimensional polyline.
Author
R. J. Purser

Definition at line 1039 of file pmat5.f90.

subroutine, public pmat5::plroti ( real(sp), dimension(3,3), intent(in)  rot3,
integer, intent(in)  n,
real(sp), dimension(n), intent(inout)  x,
real(sp), dimension(n), intent(inout)  y,
real(sp), dimension(n), intent(inout)  z 
)

Invert the rotation of a three-dimensional polyline.

Parameters
[in]rot3rotation to be inverted.
[in]nnumber of points in the polyline.
[in,out]xcartesian components of the three dimensional polyline.
[in,out]ycartesian components of the three dimensional polyline.
[in,out]zcartesian components of the three dimensional polyline.
Author
R. J. Purser

Definition at line 1061 of file pmat5.f90.

subroutine pmat5::sctoc ( real(sp), intent(in)  s,
real(sp), dimension(3), intent(inout)  xc1,
real(sp), dimension(3), intent(inout)  xc2 
)
private

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s.

Parameters
[in]sSchmidt scaling parameter
[in,out]xc1input cartesian 3-vector
[in,out]xc2output cartesian unit 3-vector
Author
R. J. Purser

Definition at line 794 of file pmat5.f90.

subroutine pmat5::sctocd ( real(sp), intent(in)  s,
real(sp), dimension(3), intent(inout)  xc1,
real(sp), dimension(3), intent(inout)  xc2,
real(sp), dimension(3,3), intent(out)  dxc2 
)

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, and its jacobian, dxc2.

Parameters
[in]sSchmidt scaling parameter.
[in,out]xc1input cartesian 3-vector.
[in,out]xc2output cartesian unit 3-vector.
[out]dxc2jacobian 1st derivative, d(xc2)/d(xc1).
Author
R. J. Purser

Definition at line 824 of file pmat5.f90.

subroutine pmat5::sctocdd ( real(sp), intent(in)  s,
real(sp), dimension(3), intent(inout)  xc1,
real(sp), dimension(3), intent(inout)  xc2,
real(sp), dimension(3,3), intent(out)  dxc2,
real(sp), dimension(3,3,3), intent(out)  ddxc2 
)

Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s, its jacobian, dxc2, and its 2nd derivative, ddxc2.

Parameters
[in]sSchmidt scaling parameter.
[in]xc1input cartesian 3-vector.
[in]xc2output cartesian unit 3-vector.
[out]dxc2jacobian 1st derivative, d(xc2)/d(xc1).
[out]ddxc22nd derivative, d^2(xc2)/(d(xc1)d(xc1)).
Author
R. J. Purser

Definition at line 865 of file pmat5.f90.

subroutine pmat5::sctog ( real(sp), dimension(3), intent(in)  xe,
real(sp), intent(out)  dlat,
real(sp), intent(out)  dlon 
)

Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Single precision version.

Parameters
[in]xeEarth-centered cartesian unit 3-vector.
[out]dlatdegrees latitude.
[out]dlondegrees longitude.
Author
R. J. Purser
Date
1994

Definition at line 372 of file pmat5.f90.

subroutine pmat5::sctogr ( real(sp), dimension(3), intent(in)  xe,
real(sp), intent(out)  rlat,
real(sp), intent(out)  rlon 
)

Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to latitude and longitude (radians) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Single precision version.

Parameters
[in]xeEarth-centered cartesian unit 3-vector.
[out]rlatradians latitude.
[out]rlonradians longitude.
Author
R. J. Purser

Definition at line 172 of file pmat5.f90.

subroutine pmat5::sframetwist ( real(sp), dimension(3), intent(in)  sxp,
real(sp), dimension(3), intent(in)  syp,
real(sp), dimension(3), intent(in)  szp,
real(sp), dimension(3), intent(in)  sxv,
real(sp), dimension(3), intent(in)  syv,
real(sp), dimension(3), intent(in)  szv,
real(sp), intent(out)  stwist 
)

Given a principal cartesian orthonormal frame, {xp,yp,zp} (i.e., at P with Earth-centered cartesians, zp), and another similar frame {xv,yv,zv} at V with Earth-centered cartesians zv, find the relative rotation angle, "twist" by which the frame at V is rotated in the counterclockwise sense relative to the parallel-transportation of P's frame to V.

Note that, by symmetry, transposing P and V leads to the opposite twist. Single precision version.

Parameters
[in]sxpP-frame cartesian X-vector.
[in]sypP-frame cartesian Y-vector.
[in]szpP-frame cartesian Z-vector.
[in]sxvV-frame cartesian X-vector.
[in]syvV-frame cartesian Y-vector.
[in]szvV-frame cartesian Z-vector.
[out]stwistrelative rotation angle (radians) of frames.
Author
R. J. Purser

Definition at line 750 of file pmat5.f90.

subroutine pmat5::sgrtoc ( real(sp), intent(in)  rlat,
real(sp), intent(in)  rlon,
real(sp), dimension(3), intent(out)  xe 
)

Transform "Geographical" to "Cartesian" coordinates.

Single proecision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
Author
R. J. Purser

Definition at line 219 of file pmat5.f90.

subroutine pmat5::sgrtocd ( real(sp), intent(in)  rlat,
real(sp), intent(in)  rlon,
real(sp), dimension(3), intent(out)  xe,
real(sp), dimension(3), intent(out)  dxedlat,
real(sp), dimension(3), intent(out)  dxedlon 
)
private

Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude.

Single precision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
[out]dxedlatDerivative, d(xe)/d(rlat).
[out]dxedlonDerivative, d(xe)/d(rlon).
Author
R. J. Purser

Definition at line 256 of file pmat5.f90.

References pmat5::grtoc::dgrtocd().

subroutine pmat5::sgrtocdd ( real(sp), intent(in)  rlat,
real(sp), intent(in)  rlon,
real(sp), dimension(3), intent(out)  xe,
real(sp), dimension(3), intent(out)  dxedlat,
real(sp), dimension(3), intent(out)  dxedlon,
real(sp), dimension(3), intent(out)  ddxedlatdlat,
real(sp), dimension(3), intent(out)  ddxedlatdlon,
real(sp), dimension(3), intent(out)  ddxedlondlon 
)

Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivatives of cartesians wrt latitude and longitude.

Single precision version.

Parameters
[in]rlatLatitude (radians) of point.
[in]rlonLongitude (radians) of point.
[out]xeEarth-centered cartesian unit 3-vector of point.
[out]dxedlatDerivative, d(xe)/d(rlat).
[out]dxedlonDerivative, d(xe)/d(rlon).
[out]ddxedlatdlatDerivative, d^2(xe)/(d(rlat)d(rlat)).
[out]ddxedlatdlonDerivative, d^2(xe)/(d(rlat)d(rlon)).
[out]ddxedlondlonDerivative, d^2(xe)/(d(rlon)d(rlon)).
Author
R. J. Purser

Definition at line 305 of file pmat5.f90.

References pmat5::gtoc::dgtocdd().

subroutine pmat5::sgtoc ( real(sp), intent(in)  dlat,
real(sp), intent(in)  dlon,
real(sp), dimension(3), intent(out)  xe 
)

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Single precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xeEarth-centered cartesian unit 3-vector.
Author
R. J. Purser
Date
1994

Definition at line 422 of file pmat5.f90.

subroutine pmat5::sgtocd ( real(sp), intent(in)  dlat,
real(sp), intent(in)  dlon,
real(sp), dimension(3), intent(out)  xe,
real(sp), dimension(3), intent(out)  dxedlat,
real(sp), dimension(3), intent(out)  dxedlon 
)

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Also, return the partial derivatives of xe wrt latitude and longitude. Single precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xecartesian unit 3-vector.
[out]dxedlatderivative, d(xe)/d(dlat).
[out]dxedlonderivative, d(xe)/d(dlon).
Author
R. J. Purser

Definition at line 469 of file pmat5.f90.

References pmat5::gtoc::dgtocd().

subroutine pmat5::sgtocdd ( real(sp), intent(in)  dlat,
real(sp), intent(in)  dlon,
real(sp), dimension(3), intent(out)  xe,
real(sp), dimension(3), intent(out)  dxedlat,
real(sp), dimension(3), intent(out)  dxedlon,
real(sp), dimension(3), intent(out)  ddxedlatdlat,
real(sp), dimension(3), intent(out)  ddxedlatdlon,
real(sp), dimension(3), intent(out)  ddxedlondlon 
)

Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to latitude and longitude (degrees) and cartesian coordinates are standard earth-centered cartesian coordinates: xe(3) pointing north, xe(1) pointing to the 0-meridian.

Also, return the 1st and 2nd partial derivatives of xe wrt latitude and longitude. Single precision version.

Parameters
[in]dlatdegrees latitude.
[in]dlondegrees longitude.
[out]xecartesian unit 3-vector.
[out]dxedlatderivative, d(xe)/d(dlat).
[out]dxedlonderivative, d(xe)/d(dlon).
[out]ddxedlatdlatderivative, d^2(xe)/(d(dlat)d(dlat)).
[out]ddxedlatdlonderivative, d^2(xe)/(d(dlat)d(dlon)).
[out]ddxedlondlonderivative, d^2(xe)/(d(dlon)d(dlon)).
Author
R. J. Purser

Definition at line 526 of file pmat5.f90.

References pmat5::gtoc::dgtocdd().

subroutine pmat5::sgtoframem ( real(sp), intent(in)  splat,
real(sp), intent(in)  splon,
real(sp), dimension(3,3), intent(out)  sorth 
)

From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location as an orthonormal matrix, orth.

Single precision version.

Parameters
[in]splatlatitude (degrees) of point.
[in]splonlongitude (degrees) of point.
[out]sorthorthonormal matrix.
Author
R. J. Purser

Definition at line 599 of file pmat5.f90.

References pmat5::gtoframe::gtoframem().

subroutine pmat5::sgtoframev ( real(sp), intent(in)  splat,
real(sp), intent(in)  splon,
real(sp), dimension(3), intent(out)  sxp,
real(sp), dimension(3), intent(out)  syp,
real(sp), dimension(3), intent(out)  szp 
)
private

Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal cartesian frame, {xp,yp,zp} in earth-centered coordinates.

Single precision version.

Parameters
[in]splatlatitude (degrees) of point.
[in]splonlongitude (degrees) of point.
[out]sxpxp unit X-basis vector of cartesian frame.
[out]sypyp unit Y-basis vector of cartesian frame.
[out]szpzp unit Z-basis vector of cartesian frame.
Author
R. J. Purser

Definition at line 636 of file pmat5.f90.

References pmat5::gtoframe::gtoframev().

subroutine pmat5::sininmap ( real(sp), intent(in)  alon0,
real(sp), intent(in)  alat0,
real(sp), dimension(3,3), intent(out)  rot3 
)
private

Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) on the projection axis.

Single precision version.

Parameters
[in]alon0geographical longitude (degrees) of projection center.
[in]alat0geographical latitude (degrees) of projection center.
[out]rot3rotation matrix.
Author
R. J. Purser
Date
1995

Definition at line 63 of file pmat5.f90.

subroutine pmat5::sinivmap ( real(sp), intent(in)  alon0,
real(sp), intent(in)  alat0,
real(sp), dimension(3,3), intent(out)  rot3 
)

Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components to the alternative cartesian frame oriented so as to put geographical point (ALAT0,ALON0) at the viewing nadir.

Single precision version.

Parameters
[in]alon0geographical longitude (degrees) of viewing nadir.
[in]alat0geographical latitude (degrees) of viewing nadir.
[out]rot3rotation matrix.
Author
R. J. Purser
Date
1995

Definition at line 107 of file pmat5.f90.

subroutine pmat5::sparaframe ( real(sp), dimension(3), intent(in)  sxp,
real(sp), dimension(3), intent(in)  syp,
real(sp), dimension(3), intent(in)  szp,
real(sp), dimension(3), intent(out)  sxv,
real(sp), dimension(3), intent(out)  syv,
real(sp), dimension(3), intent(in)  szv 
)

Take a principal reference orthonormal frame, {xp,yp,zp} and a dependent point defined by unit vector, zv, and complete the V-frame cartesian components, {xv,yv}, that are the result of parallel-transport of {xp,yp} along the geodesic between P and V.

Single precision version.

Parameters
[in]sxpreference orthonormal P-frame cartesian X-vector.
[in]sypreference orthonormal P-frame cartesian Y-vector.
[in]szpreference orthonormal P-frame cartesian Z-vector.
[out]sxvV-frame cartesian X-vector.
[out]syvV-frame cartesian Y-vector.
[in]szvdependent point zenith, V-frame cartesian Z-vector.
Author
R. J. Purser

Definition at line 695 of file pmat5.f90.


The documentation for this module was generated from the following file: