/src/proj/src/conversions/cart.cpp
Line | Count | Source (jump to first uncovered line) |
1 | | /****************************************************************************** |
2 | | * Project: PROJ.4 |
3 | | * Purpose: Convert between ellipsoidal, geodetic coordinates and |
4 | | * cartesian, geocentric coordinates. |
5 | | * |
6 | | * Formally, this functionality is also found in the PJ_geocent.c |
7 | | * code. |
8 | | * |
9 | | * Actually, however, the PJ_geocent transformations are carried |
10 | | * out in concert between 2D stubs in PJ_geocent.c and 3D code |
11 | | * placed in pj_transform.c. |
12 | | * |
13 | | * For pipeline-style datum shifts, we do need direct access |
14 | | * to the full 3D interface for this functionality. |
15 | | * |
16 | | * Hence this code, which may look like "just another PJ_geocent" |
17 | | * but really is something substantially different. |
18 | | * |
19 | | * Author: Thomas Knudsen, thokn@sdfe.dk |
20 | | * |
21 | | ****************************************************************************** |
22 | | * Copyright (c) 2016, Thomas Knudsen / SDFE |
23 | | * |
24 | | * Permission is hereby granted, free of charge, to any person obtaining a |
25 | | * copy of this software and associated documentation files (the "Software"), |
26 | | * to deal in the Software without restriction, including without limitation |
27 | | * the rights to use, copy, modify, merge, publish, distribute, sublicense, |
28 | | * and/or sell copies of the Software, and to permit persons to whom the |
29 | | * Software is furnished to do so, subject to the following conditions: |
30 | | * |
31 | | * The above copyright notice and this permission notice shall be included |
32 | | * in all copies or substantial portions of the Software. |
33 | | * |
34 | | * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS |
35 | | * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, |
36 | | * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL |
37 | | * THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER |
38 | | * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING |
39 | | * FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER |
40 | | * DEALINGS IN THE SOFTWARE. |
41 | | *****************************************************************************/ |
42 | | |
43 | | #include "proj_internal.h" |
44 | | #include <math.h> |
45 | | |
46 | | PROJ_HEAD(cart, "Geodetic/cartesian conversions"); |
47 | | |
48 | | /************************************************************** |
49 | | CARTESIAN / GEODETIC CONVERSIONS |
50 | | *************************************************************** |
51 | | This material follows: |
52 | | |
53 | | Bernhard Hofmann-Wellenhof & Helmut Moritz: |
54 | | Physical Geodesy, 2nd edition. |
55 | | Springer, 2005. |
56 | | |
57 | | chapter 5.6: Coordinate transformations |
58 | | (HM, below), |
59 | | |
60 | | and |
61 | | |
62 | | Wikipedia: Geographic Coordinate Conversion, |
63 | | https://en.wikipedia.org/wiki/Geographic_coordinate_conversion |
64 | | |
65 | | (WP, below). |
66 | | |
67 | | The cartesian-to-geodetic conversion is based on Bowring's |
68 | | celebrated method: |
69 | | |
70 | | B. R. Bowring: |
71 | | Transformation from spatial to geographical coordinates |
72 | | Survey Review 23(181), pp. 323-327, 1976 |
73 | | |
74 | | (BB, below), |
75 | | |
76 | | but could probably use some TLC from a newer and faster |
77 | | algorithm: |
78 | | |
79 | | Toshio Fukushima: |
80 | | Transformation from Cartesian to Geodetic Coordinates |
81 | | Accelerated by Halley’s Method |
82 | | Journal of Geodesy, February 2006 |
83 | | |
84 | | (TF, below). |
85 | | |
86 | | Close to the poles, we avoid singularities by switching to an |
87 | | approximation requiring knowledge of the geocentric radius |
88 | | at the given latitude. For this, we use an adaptation of the |
89 | | formula given in: |
90 | | |
91 | | Wikipedia: Earth Radius |
92 | | https://en.wikipedia.org/wiki/Earth_radius#Radius_at_a_given_geodetic_latitude |
93 | | (Derivation and commentary at https://gis.stackexchange.com/q/20200) |
94 | | |
95 | | (WP2, below) |
96 | | |
97 | | These routines are probably not as robust at those in |
98 | | geocent.c, at least they haven't been through as heavy |
99 | | use as their geocent sisters. Some care has been taken |
100 | | to avoid singularities, but extreme cases (e.g. setting |
101 | | es, the squared eccentricity, to 1), will cause havoc. |
102 | | |
103 | | **************************************************************/ |
104 | | |
105 | | /*********************************************************************/ |
106 | 277 | static double normal_radius_of_curvature(double a, double es, double sinphi) { |
107 | | /*********************************************************************/ |
108 | 277 | if (es == 0) |
109 | 140 | return a; |
110 | | /* This is from WP. HM formula 2-149 gives an a,b version */ |
111 | 137 | return a / sqrt(1 - es * sinphi * sinphi); |
112 | 277 | } |
113 | | |
114 | | /*********************************************************************/ |
115 | | static double geocentric_radius(double a, double b_div_a, double cosphi, |
116 | 19 | double sinphi) { |
117 | | /********************************************************************* |
118 | | Return the geocentric radius at latitude phi, of an ellipsoid |
119 | | with semimajor axis a and semiminor axis b. |
120 | | |
121 | | This is from WP2, but uses hypot() for potentially better |
122 | | numerical robustness |
123 | | ***********************************************************************/ |
124 | | // Non-optimized version: |
125 | | // const double b = a * b_div_a; |
126 | | // return hypot(a * a * cosphi, b * b * sinphi) / |
127 | | // hypot(a * cosphi, b * sinphi); |
128 | 19 | const double cosphi_squared = cosphi * cosphi; |
129 | 19 | const double sinphi_squared = sinphi * sinphi; |
130 | 19 | const double b_div_a_squared = b_div_a * b_div_a; |
131 | 19 | const double b_div_a_squared_mul_sinphi_squared = |
132 | 19 | b_div_a_squared * sinphi_squared; |
133 | 19 | return a * sqrt((cosphi_squared + |
134 | 19 | b_div_a_squared * b_div_a_squared_mul_sinphi_squared) / |
135 | 19 | (cosphi_squared + b_div_a_squared_mul_sinphi_squared)); |
136 | 19 | } |
137 | | |
138 | | /*********************************************************************/ |
139 | 277 | static PJ_XYZ cartesian(PJ_LPZ geod, PJ *P) { |
140 | | /*********************************************************************/ |
141 | 277 | PJ_XYZ xyz; |
142 | | |
143 | 277 | const double cosphi = cos(geod.phi); |
144 | 277 | const double sinphi = sin(geod.phi); |
145 | 277 | const double N = normal_radius_of_curvature(P->a, P->es, sinphi); |
146 | | |
147 | | /* HM formula 5-27 (z formula follows WP) */ |
148 | 277 | xyz.x = (N + geod.z) * cosphi * cos(geod.lam); |
149 | 277 | xyz.y = (N + geod.z) * cosphi * sin(geod.lam); |
150 | 277 | xyz.z = (N * (1 - P->es) + geod.z) * sinphi; |
151 | | |
152 | 277 | return xyz; |
153 | 277 | } |
154 | | |
155 | | /*********************************************************************/ |
156 | 19 | static PJ_LPZ geodetic(PJ_XYZ cart, PJ *P) { |
157 | | /*********************************************************************/ |
158 | 19 | PJ_LPZ lpz; |
159 | | |
160 | | // Normalize (x,y,z) to the unit sphere/ellipsoid. |
161 | | #if (defined(__i386__) && !defined(__SSE__)) || defined(_M_IX86) |
162 | | // i386 (actually non-SSE) code path to make following test case of |
163 | | // testvarious happy |
164 | | // "echo 6378137.00 -0.00 0.00 | bin/cs2cs +proj=geocent +datum=WGS84 +to |
165 | | // +proj=latlong +datum=WGS84" |
166 | | const double x_div_a = cart.x / P->a; |
167 | | const double y_div_a = cart.y / P->a; |
168 | | const double z_div_a = cart.z / P->a; |
169 | | #else |
170 | 19 | const double x_div_a = cart.x * P->ra; |
171 | 19 | const double y_div_a = cart.y * P->ra; |
172 | 19 | const double z_div_a = cart.z * P->ra; |
173 | 19 | #endif |
174 | | |
175 | | /* Perpendicular distance from point to Z-axis (HM eq. 5-28) */ |
176 | 19 | const double p_div_a = sqrt(x_div_a * x_div_a + y_div_a * y_div_a); |
177 | | |
178 | | #if 0 |
179 | | /* HM eq. (5-37) */ |
180 | | const double theta = atan2 (cart.z * P->a, p * P->b); |
181 | | |
182 | | /* HM eq. (5-36) (from BB, 1976) */ |
183 | | const double c = cos(theta); |
184 | | const double s = sin(theta); |
185 | | #else |
186 | 19 | const double b_div_a = 1 - P->f; // = P->b / P->a |
187 | 19 | const double p_div_a_b_div_a = p_div_a * b_div_a; |
188 | 19 | const double norm = |
189 | 19 | sqrt(z_div_a * z_div_a + p_div_a_b_div_a * p_div_a_b_div_a); |
190 | 19 | double c, s; |
191 | 19 | if (norm != 0) { |
192 | 0 | const double inv_norm = 1.0 / norm; |
193 | 0 | c = p_div_a_b_div_a * inv_norm; |
194 | 0 | s = z_div_a * inv_norm; |
195 | 19 | } else { |
196 | 19 | c = 1; |
197 | 19 | s = 0; |
198 | 19 | } |
199 | 19 | #endif |
200 | | |
201 | 19 | const double y_phi = z_div_a + P->e2s * b_div_a * s * s * s; |
202 | 19 | const double x_phi = p_div_a - P->es * c * c * c; |
203 | 19 | const double norm_phi = sqrt(y_phi * y_phi + x_phi * x_phi); |
204 | 19 | double cosphi, sinphi; |
205 | 19 | if (norm_phi != 0) { |
206 | 19 | const double inv_norm_phi = 1.0 / norm_phi; |
207 | 19 | cosphi = x_phi * inv_norm_phi; |
208 | 19 | sinphi = y_phi * inv_norm_phi; |
209 | 19 | } else { |
210 | 0 | cosphi = 1; |
211 | 0 | sinphi = 0; |
212 | 0 | } |
213 | 19 | if (x_phi <= 0) { |
214 | | // this happen on non-sphere ellipsoid when x,y,z is very close to 0 |
215 | | // there is no single solution to the cart->geodetic conversion in |
216 | | // that case, clamp to -90/90 deg and avoid a discontinuous boundary |
217 | | // near the poles |
218 | 19 | lpz.phi = cart.z >= 0 ? M_HALFPI : -M_HALFPI; |
219 | 19 | cosphi = 0; |
220 | 19 | sinphi = cart.z >= 0 ? 1 : -1; |
221 | 19 | } else { |
222 | 0 | lpz.phi = atan(y_phi / x_phi); |
223 | 0 | } |
224 | 19 | lpz.lam = atan2(y_div_a, x_div_a); |
225 | | |
226 | 19 | if (cosphi < 1e-6) { |
227 | | /* poleward of 89.99994 deg, we avoid division by zero */ |
228 | | /* by computing the height as the cartesian z value */ |
229 | | /* minus the geocentric radius of the Earth at the given */ |
230 | | /* latitude */ |
231 | 19 | const double r = geocentric_radius(P->a, b_div_a, cosphi, sinphi); |
232 | 19 | lpz.z = fabs(cart.z) - r; |
233 | 19 | } else { |
234 | 0 | const double N = normal_radius_of_curvature(P->a, P->es, sinphi); |
235 | 0 | lpz.z = P->a * p_div_a / cosphi - N; |
236 | 0 | } |
237 | | |
238 | 19 | return lpz; |
239 | 19 | } |
240 | | |
241 | | /* In effect, 2 cartesian coordinates of a point on the ellipsoid. Rather |
242 | | * pointless, but... */ |
243 | 0 | static PJ_XY cart_forward(PJ_LP lp, PJ *P) { |
244 | 0 | PJ_COORD point; |
245 | 0 | point.lp = lp; |
246 | 0 | point.lpz.z = 0; |
247 | |
|
248 | 0 | const auto xyz = cartesian(point.lpz, P); |
249 | 0 | point.xyz = xyz; |
250 | 0 | return point.xy; |
251 | 0 | } |
252 | | |
253 | | /* And the other way round. Still rather pointless, but... */ |
254 | 0 | static PJ_LP cart_reverse(PJ_XY xy, PJ *P) { |
255 | 0 | PJ_COORD point; |
256 | 0 | point.xy = xy; |
257 | 0 | point.xyz.z = 0; |
258 | |
|
259 | 0 | const auto lpz = geodetic(point.xyz, P); |
260 | 0 | point.lpz = lpz; |
261 | 0 | return point.lp; |
262 | 0 | } |
263 | | |
264 | | /*********************************************************************/ |
265 | 19.1k | PJ *PJ_CONVERSION(cart, 1) { |
266 | | /*********************************************************************/ |
267 | 19.1k | P->fwd3d = cartesian; |
268 | 19.1k | P->inv3d = geodetic; |
269 | 19.1k | P->fwd = cart_forward; |
270 | 19.1k | P->inv = cart_reverse; |
271 | 19.1k | P->left = PJ_IO_UNITS_RADIANS; |
272 | 19.1k | P->right = PJ_IO_UNITS_CARTESIAN; |
273 | 19.1k | return P; |
274 | 19.1k | } |