grid_tools
1.9.0
|
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 |
Functions/Subroutines | |
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 | 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 | 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 | 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 | 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 | 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, 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 | 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, 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 | sctoc (s, xc1, xc2) |
Evaluate Schmidt transformation, xc1 –> xc2, with scaling parameter s. 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 | sgrtocd (rlat, rlon, xe, dxedlat, dxedlon) |
Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude. 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 | 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... | |
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... | |
Utility routines for orienting the globe and basic geographical mappings.
subroutine pmat5::dctoc | ( | real(dp), intent(in) | s, |
real(dp), dimension(3), intent(inout) | xc1, | ||
real(dp), dimension(3), intent(inout) | xc2 | ||
) |
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.
[in] | s | Schmidt scaling parameter. |
[in,out] | xc1 | input cartesian 3-vector. |
[in,out] | xc2 | output cartesian unit 3-vector. |
[out] | dxc2 | jacobian 1st derivative, d(xc2)/d(xc1). |
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.
[in] | s | Schmidt scaling parameter. |
[in,out] | xc1 | input cartesian 3-vector. |
[in,out] | xc2 | output cartesian unit 3-vector. |
[out] | dxc2 | jacobian 1st derivative, d(xc2)/d(xc1). |
[out] | ddxc2 | 2nd derivative, d^2(xc2)/(d(xc1)d(xc1)). |
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.
[in] | xe | Earth-centered cartesian unit 3-vector. |
[out] | dlat | degrees latitude. |
[out] | dlon | degrees longitude. |
Definition at line 398 of file pmat5.f90.
References pietc::rtod, and pietc::u0.
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.
[in] | xe | Earth-centered cartesian unit 3-vector. |
[out] | rlat | radians latitude. |
[out] | rlon | radians longitude. |
Definition at line 198 of file pmat5.f90.
References pietc::u0.
|
private |
|
private |
Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude.
Double precision version.
[in] | rlat | Latitude (radians) of point. |
[in] | rlon | Longitude (radians) of point. |
[out] | xe | Earth-centered cartesian unit 3-vector of point. |
[out] | dxedlat | Derivative, d(xe)/d(rlat). |
[out] | dxedlon | Derivative, d(xe)/d(rlon). |
Definition at line 280 of file pmat5.f90.
References pietc::u0.
Referenced by sgrtocd().
|
private |
Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivatives of cartesians wrt latitude and longitude.
Double precision version.
[in] | rlat | Latitude (radians) of point. |
[in] | rlon | Longitude (radians) of point. |
[out] | xe | Earth-centered cartesian unit 3-vector of point. |
[out] | dxedlat | Derivative, d(xe)/d(rlat). |
[out] | dxedlon | Derivative, d(xe)/d(rlon). |
[out] | ddxedlatdlat | Derivative, d^2(xe)/(d(rlat)d(rlat)). |
[out] | ddxedlatdlon | Derivative, d^2(xe)/(d(rlat)d(rlon)). |
[out] | ddxedlondlon | Derivative, d^2(xe)/(d(rlon)d(rlon)). |
Definition at line 340 of file pmat5.f90.
References pietc::u0.
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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | cartesian unit 3-vector. |
Definition at line 445 of file pmat5.f90.
References pietc::dtor.
|
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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | cartesian unit 3-vector. |
[out] | dxedlat | derivative, d(xe)/d(dlat). |
[out] | dxedlon | derivative, d(xe)/d(dlon). |
Definition at line 496 of file pmat5.f90.
References pietc::dtor, and pietc::u0.
Referenced by sgtocd().
|
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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | cartesian unit 3-vector. |
[out] | dxedlat | d(xe)/d(dlat). |
[out] | dxedlon | d(xe)/d(dlon). |
[out] | ddxedlatdlat | derivative, d^2(xe)/(d(dlat)d(dlat)). |
[out] | ddxedlatdlon | derivative, d^2(xe)/(d(dlat)d(dlon)). |
[out] | ddxedlondlon | derivative, d^2(xe)/(d(dlon)d(dlon)). |
Definition at line 565 of file pmat5.f90.
References pietc::dtor, and pietc::u0.
Referenced by sgrtocdd(), and sgtocdd().
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.
[in] | alon0 | geographical longitude (degrees) of projection center. |
[in] | alat0 | geographical latitude (degrees) of projection center. |
[out] | rot3 | rotation matrix. |
Definition at line 86 of file pmat5.f90.
References pietc::dtor, and pietc::u0.
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.
[in] | alon0 | geographical longitude (degrees) of viewing nadir. |
[in] | alat0 | geographical latitude (degrees) of viewing nadir. |
[out] | rot3 | rotation matrix. |
Definition at line 140 of file pmat5.f90.
References pietc::dtor, and pietc::u0.
|
private |
Apply a constant rotation to a three dimensional polyline.
[in] | rot3 | rotation matrix. |
[in] | n | number of points in the polyline. |
[in,out] | x | cartesian components of the three dimensional polyline. |
[in,out] | y | cartesian components of the three dimensional polyline. |
[in,out] | z | cartesian components of the three dimensional polyline. |
|
private |
Invert the rotation of a three-dimensional polyline.
[in] | rot3 | rotation to be inverted. |
[in] | n | number of points in the polyline. |
[in,out] | x | cartesian components of the three dimensional polyline. |
[in,out] | y | cartesian components of the three dimensional polyline. |
[in,out] | z | cartesian components of the three dimensional polyline. |
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.
[in] | xp | P-frame cartesian X-vector. |
[in] | yp | P-frame cartesian Y-vector. |
[in] | zp | P-frame cartesian Z-vector. |
[in] | xv | V-frame cartesian X-vector. |
[in] | yv | V-frame cartesian Y-vector. |
[in] | zv | V-frame cartesian Z-vector. |
[out] | twist | relative rotation angle (radians) of frames. |
|
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.
[in] | plat | latitude (degrees) of point. |
[in] | plon | longitude (degrees) of point. |
[out] | orth | orthonormal matrix. |
Definition at line 617 of file pmat5.f90.
References gtoframev().
Referenced by sgtoframem().
|
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.
[in] | plat | latitude (degrees) of point. |
[in] | plon | longitude (degrees) of point. |
[out] | xp | unit X-basis vector of cartesian frame. |
[out] | yp | unit Y-basis vector of cartesian frame. |
[out] | zp | unit Z-basis vector of cartesian frame. |
Definition at line 662 of file pmat5.f90.
References pietc::u0, and pietc::u1.
Referenced by gtoframem(), and sgtoframev().
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.
[in] | xp | reference orthonormal P-frame cartesian X-vector. |
[in] | yp | reference orthonormal P-frame cartesian Y-vector. |
[in] | zp | reference orthonormal P-frame cartesian Z-vector. |
[out] | xv | V-frame cartesian X-vector. |
[out] | yv | V-frame cartesian Y-vector. |
[in] | zv | dependent point zenith, V-frame cartesian Z-vector. |
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.
[in] | s | Schmidt scaling parameter. |
[in] | n | number of points in the polyline. |
[in,out] | x | cartesian components of the three dimensional polyline. |
[in,out] | y | cartesian components of the three dimensional polyline. |
[in,out] | z | cartesian components of the three dimensional polyline. |
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.
[in] | rot3 | rotation matrix. |
[in] | n | number of points in the polyline. |
[in,out] | x | cartesian components of the three dimensional polyline. |
[in,out] | y | cartesian components of the three dimensional polyline. |
[in,out] | z | cartesian components of the three dimensional polyline. |
Definition at line 1040 of file pmat5.f90.
References pietc::t.
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.
[in] | rot3 | rotation to be inverted. |
[in] | n | number of points in the polyline. |
[in,out] | x | cartesian components of the three dimensional polyline. |
[in,out] | y | cartesian components of the three dimensional polyline. |
[in,out] | z | cartesian components of the three dimensional polyline. |
Definition at line 1062 of file pmat5.f90.
References pietc::t.
|
private |
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.
[in] | s | Schmidt scaling parameter. |
[in,out] | xc1 | input cartesian 3-vector. |
[in,out] | xc2 | output cartesian unit 3-vector. |
[out] | dxc2 | jacobian 1st derivative, d(xc2)/d(xc1). |
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.
[in] | s | Schmidt scaling parameter. |
[in] | xc1 | input cartesian 3-vector. |
[in] | xc2 | output cartesian unit 3-vector. |
[out] | dxc2 | jacobian 1st derivative, d(xc2)/d(xc1). |
[out] | ddxc2 | 2nd derivative, d^2(xc2)/(d(xc1)d(xc1)). |
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.
[in] | xe | Earth-centered cartesian unit 3-vector. |
[out] | dlat | degrees latitude. |
[out] | dlon | degrees longitude. |
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.
[in] | xe | Earth-centered cartesian unit 3-vector. |
[out] | rlat | radians latitude. |
[out] | rlon | radians longitude. |
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.
[in] | sxp | P-frame cartesian X-vector. |
[in] | syp | P-frame cartesian Y-vector. |
[in] | szp | P-frame cartesian Z-vector. |
[in] | sxv | V-frame cartesian X-vector. |
[in] | syv | V-frame cartesian Y-vector. |
[in] | szv | V-frame cartesian Z-vector. |
[out] | stwist | relative rotation angle (radians) of frames. |
subroutine pmat5::sgrtoc | ( | real(sp), intent(in) | rlat, |
real(sp), intent(in) | rlon, | ||
real(sp), dimension(3), intent(out) | xe | ||
) |
|
private |
Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesians wrt latitude and longitude.
Single precision version.
[in] | rlat | Latitude (radians) of point. |
[in] | rlon | Longitude (radians) of point. |
[out] | xe | Earth-centered cartesian unit 3-vector of point. |
[out] | dxedlat | Derivative, d(xe)/d(rlat). |
[out] | dxedlon | Derivative, d(xe)/d(rlon). |
Definition at line 257 of file pmat5.f90.
References 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.
[in] | rlat | Latitude (radians) of point. |
[in] | rlon | Longitude (radians) of point. |
[out] | xe | Earth-centered cartesian unit 3-vector of point. |
[out] | dxedlat | Derivative, d(xe)/d(rlat). |
[out] | dxedlon | Derivative, d(xe)/d(rlon). |
[out] | ddxedlatdlat | Derivative, d^2(xe)/(d(rlat)d(rlat)). |
[out] | ddxedlatdlon | Derivative, d^2(xe)/(d(rlat)d(rlon)). |
[out] | ddxedlondlon | Derivative, d^2(xe)/(d(rlon)d(rlon)). |
Definition at line 307 of file pmat5.f90.
References 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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | Earth-centered cartesian unit 3-vector. |
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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | cartesian unit 3-vector. |
[out] | dxedlat | derivative, d(xe)/d(dlat). |
[out] | dxedlon | derivative, d(xe)/d(dlon). |
Definition at line 470 of file pmat5.f90.
References 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.
[in] | dlat | degrees latitude. |
[in] | dlon | degrees longitude. |
[out] | xe | cartesian unit 3-vector. |
[out] | dxedlat | derivative, d(xe)/d(dlat). |
[out] | dxedlon | derivative, d(xe)/d(dlon). |
[out] | ddxedlatdlat | derivative, d^2(xe)/(d(dlat)d(dlat)). |
[out] | ddxedlatdlon | derivative, d^2(xe)/(d(dlat)d(dlon)). |
[out] | ddxedlondlon | derivative, d^2(xe)/(d(dlon)d(dlon)). |
Definition at line 528 of file pmat5.f90.
References 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.
[in] | splat | latitude (degrees) of point. |
[in] | splon | longitude (degrees) of point. |
[out] | sorth | orthonormal matrix. |
Definition at line 600 of file pmat5.f90.
References gtoframem().
|
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.
[in] | splat | latitude (degrees) of point. |
[in] | splon | longitude (degrees) of point. |
[out] | sxp | xp unit X-basis vector of cartesian frame. |
[out] | syp | yp unit Y-basis vector of cartesian frame. |
[out] | szp | zp unit Z-basis vector of cartesian frame. |
Definition at line 637 of file pmat5.f90.
References gtoframev().
|
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.
[in] | alon0 | geographical longitude (degrees) of projection center. |
[in] | alat0 | geographical latitude (degrees) of projection center. |
[out] | rot3 | rotation matrix. |
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.
[in] | alon0 | geographical longitude (degrees) of viewing nadir. |
[in] | alat0 | geographical latitude (degrees) of viewing nadir. |
[out] | rot3 | rotation matrix. |
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.
[in] | sxp | reference orthonormal P-frame cartesian X-vector. |
[in] | syp | reference orthonormal P-frame cartesian Y-vector. |
[in] | szp | reference orthonormal P-frame cartesian Z-vector. |
[out] | sxv | V-frame cartesian X-vector. |
[out] | syv | V-frame cartesian Y-vector. |
[in] | szv | dependent point zenith, V-frame cartesian Z-vector. |