/src/proj/src/projections/omerc.cpp
Line | Count | Source (jump to first uncovered line) |
1 | | /* |
2 | | ** Copyright (c) 2003, 2006 Gerald I. Evenden |
3 | | */ |
4 | | /* |
5 | | ** Permission is hereby granted, free of charge, to any person obtaining |
6 | | ** a copy of this software and associated documentation files (the |
7 | | ** "Software"), to deal in the Software without restriction, including |
8 | | ** without limitation the rights to use, copy, modify, merge, publish, |
9 | | ** distribute, sublicense, and/or sell copies of the Software, and to |
10 | | ** permit persons to whom the Software is furnished to do so, subject to |
11 | | ** the following conditions: |
12 | | ** |
13 | | ** The above copyright notice and this permission notice shall be |
14 | | ** included in all copies or substantial portions of the Software. |
15 | | ** |
16 | | ** THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, |
17 | | ** EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF |
18 | | ** MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. |
19 | | ** IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY |
20 | | ** CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, |
21 | | ** TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE |
22 | | ** SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE. |
23 | | */ |
24 | | |
25 | | #include <errno.h> |
26 | | #include <math.h> |
27 | | |
28 | | #include "proj.h" |
29 | | #include "proj_internal.h" |
30 | | |
31 | | PROJ_HEAD(omerc, "Oblique Mercator") |
32 | | "\n\tCyl, Sph&Ell no_rot\n\t" |
33 | | "alpha= [gamma=] [no_off] lonc= or\n\t lon_1= lat_1= lon_2= lat_2="; |
34 | | |
35 | | namespace { // anonymous namespace |
36 | | struct pj_omerc_data { |
37 | | double A, B, E, AB, ArB, BrA, rB, singam, cosgam, sinrot, cosrot; |
38 | | double v_pole_n, v_pole_s, u_0; |
39 | | int no_rot; |
40 | | }; |
41 | | } // anonymous namespace |
42 | | |
43 | 2.42k | #define TOL 1.e-7 |
44 | 824 | #define EPS 1.e-10 |
45 | | |
46 | 0 | static PJ_XY omerc_e_forward(PJ_LP lp, PJ *P) { /* Ellipsoidal, forward */ |
47 | 0 | PJ_XY xy = {0.0, 0.0}; |
48 | 0 | struct pj_omerc_data *Q = static_cast<struct pj_omerc_data *>(P->opaque); |
49 | 0 | double u, v; |
50 | |
|
51 | 0 | if (fabs(fabs(lp.phi) - M_HALFPI) > EPS) { |
52 | 0 | const double W = Q->E / pow(pj_tsfn(lp.phi, sin(lp.phi), P->e), Q->B); |
53 | 0 | const double one_div_W = 1. / W; |
54 | 0 | const double S = .5 * (W - one_div_W); |
55 | 0 | const double T = .5 * (W + one_div_W); |
56 | 0 | const double V = sin(Q->B * lp.lam); |
57 | 0 | const double U = (S * Q->singam - V * Q->cosgam) / T; |
58 | 0 | if (fabs(fabs(U) - 1.0) < EPS) { |
59 | 0 | proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN); |
60 | 0 | return xy; |
61 | 0 | } |
62 | 0 | v = 0.5 * Q->ArB * log((1. - U) / (1. + U)); |
63 | 0 | const double temp = cos(Q->B * lp.lam); |
64 | 0 | if (fabs(temp) < TOL) { |
65 | 0 | u = Q->A * lp.lam; |
66 | 0 | } else { |
67 | 0 | u = Q->ArB * atan2((S * Q->cosgam + V * Q->singam), temp); |
68 | 0 | } |
69 | 0 | } else { |
70 | 0 | v = lp.phi > 0 ? Q->v_pole_n : Q->v_pole_s; |
71 | 0 | u = Q->ArB * lp.phi; |
72 | 0 | } |
73 | 0 | if (Q->no_rot) { |
74 | 0 | xy.x = u; |
75 | 0 | xy.y = v; |
76 | 0 | } else { |
77 | 0 | u -= Q->u_0; |
78 | 0 | xy.x = v * Q->cosrot + u * Q->sinrot; |
79 | 0 | xy.y = u * Q->cosrot - v * Q->sinrot; |
80 | 0 | } |
81 | 0 | return xy; |
82 | 0 | } |
83 | | |
84 | 0 | static PJ_LP omerc_e_inverse(PJ_XY xy, PJ *P) { /* Ellipsoidal, inverse */ |
85 | 0 | PJ_LP lp = {0.0, 0.0}; |
86 | 0 | struct pj_omerc_data *Q = static_cast<struct pj_omerc_data *>(P->opaque); |
87 | 0 | double u, v, Qp, Sp, Tp, Vp, Up; |
88 | |
|
89 | 0 | if (Q->no_rot) { |
90 | 0 | v = xy.y; |
91 | 0 | u = xy.x; |
92 | 0 | } else { |
93 | 0 | v = xy.x * Q->cosrot - xy.y * Q->sinrot; |
94 | 0 | u = xy.y * Q->cosrot + xy.x * Q->sinrot + Q->u_0; |
95 | 0 | } |
96 | 0 | Qp = exp(-Q->BrA * v); |
97 | 0 | if (Qp == 0) { |
98 | 0 | proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN); |
99 | 0 | return proj_coord_error().lp; |
100 | 0 | } |
101 | 0 | Sp = .5 * (Qp - 1. / Qp); |
102 | 0 | Tp = .5 * (Qp + 1. / Qp); |
103 | 0 | Vp = sin(Q->BrA * u); |
104 | 0 | Up = (Vp * Q->cosgam + Sp * Q->singam) / Tp; |
105 | 0 | if (fabs(fabs(Up) - 1.) < EPS) { |
106 | 0 | lp.lam = 0.; |
107 | 0 | lp.phi = Up < 0. ? -M_HALFPI : M_HALFPI; |
108 | 0 | } else { |
109 | 0 | lp.phi = Q->E / sqrt((1. + Up) / (1. - Up)); |
110 | 0 | if ((lp.phi = pj_phi2(P->ctx, pow(lp.phi, 1. / Q->B), P->e)) == |
111 | 0 | HUGE_VAL) { |
112 | 0 | proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN); |
113 | 0 | return lp; |
114 | 0 | } |
115 | 0 | lp.lam = |
116 | 0 | -Q->rB * atan2((Sp * Q->cosgam - Vp * Q->singam), cos(Q->BrA * u)); |
117 | 0 | } |
118 | 0 | return lp; |
119 | 0 | } |
120 | | |
121 | 842 | PJ *PJ_PROJECTION(omerc) { |
122 | 842 | double con, com, cosph0, D, F, H, L, sinph0, p, J, |
123 | 842 | gamma = 0, gamma0, lamc = 0, lam1 = 0, lam2 = 0, phi1 = 0, phi2 = 0, |
124 | 842 | alpha_c = 0; |
125 | 842 | int alp, gam, no_off = 0; |
126 | | |
127 | 842 | struct pj_omerc_data *Q = static_cast<struct pj_omerc_data *>( |
128 | 842 | calloc(1, sizeof(struct pj_omerc_data))); |
129 | 842 | if (nullptr == Q) |
130 | 0 | return pj_default_destructor(P, PROJ_ERR_OTHER /*ENOMEM*/); |
131 | 842 | P->opaque = Q; |
132 | | |
133 | 842 | Q->no_rot = pj_param(P->ctx, P->params, "bno_rot").i; |
134 | 842 | if ((alp = pj_param(P->ctx, P->params, "talpha").i) != 0) |
135 | 384 | alpha_c = pj_param(P->ctx, P->params, "ralpha").f; |
136 | 842 | if ((gam = pj_param(P->ctx, P->params, "tgamma").i) != 0) |
137 | 154 | gamma = pj_param(P->ctx, P->params, "rgamma").f; |
138 | 842 | if (alp || gam) { |
139 | 437 | lamc = pj_param(P->ctx, P->params, "rlonc").f; |
140 | 437 | no_off = |
141 | | /* For libproj4 compatibility */ |
142 | 437 | pj_param(P->ctx, P->params, "tno_off").i |
143 | | /* for backward compatibility */ |
144 | 437 | || pj_param(P->ctx, P->params, "tno_uoff").i; |
145 | 437 | if (no_off) { |
146 | | /* Mark the parameter as used, so that the pj_get_def() return them |
147 | | */ |
148 | 74 | pj_param(P->ctx, P->params, "sno_uoff"); |
149 | 74 | pj_param(P->ctx, P->params, "sno_off"); |
150 | 74 | } |
151 | 437 | } else { |
152 | 405 | lam1 = pj_param(P->ctx, P->params, "rlon_1").f; |
153 | 405 | phi1 = pj_param(P->ctx, P->params, "rlat_1").f; |
154 | 405 | lam2 = pj_param(P->ctx, P->params, "rlon_2").f; |
155 | 405 | phi2 = pj_param(P->ctx, P->params, "rlat_2").f; |
156 | 405 | con = fabs(phi1); |
157 | | |
158 | 405 | if (fabs(phi1) > M_HALFPI - TOL) { |
159 | 3 | proj_log_error( |
160 | 3 | P, _("Invalid value for lat_1: |lat_1| should be < 90°")); |
161 | 3 | return pj_default_destructor(P, |
162 | 3 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
163 | 3 | } |
164 | 402 | if (fabs(phi2) > M_HALFPI - TOL) { |
165 | 3 | proj_log_error( |
166 | 3 | P, _("Invalid value for lat_2: |lat_2| should be < 90°")); |
167 | 3 | return pj_default_destructor(P, |
168 | 3 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
169 | 3 | } |
170 | | |
171 | 399 | if (fabs(phi1 - phi2) <= TOL) { |
172 | 7 | proj_log_error(P, |
173 | 7 | _("Invalid value for lat_1/lat_2: lat_1 should be " |
174 | 7 | "different from lat_2")); |
175 | 7 | return pj_default_destructor(P, |
176 | 7 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
177 | 7 | } |
178 | | |
179 | 392 | if (con <= TOL) { |
180 | 4 | proj_log_error( |
181 | 4 | P, |
182 | 4 | _("Invalid value for lat_1: lat_1 should be different from 0")); |
183 | 4 | return pj_default_destructor(P, |
184 | 4 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
185 | 4 | } |
186 | | |
187 | 388 | if (fabs(fabs(P->phi0) - M_HALFPI) <= TOL) { |
188 | 1 | proj_log_error( |
189 | 1 | P, _("Invalid value for lat_0: |lat_0| should be < 90°")); |
190 | 1 | return pj_default_destructor(P, |
191 | 1 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
192 | 1 | } |
193 | 388 | } |
194 | | |
195 | 824 | if (pj_param(P->ctx, P->params, "rlon_0").i) { |
196 | 122 | proj_log_trace(P, _("lon_0 is ignored.")); |
197 | 122 | } |
198 | | |
199 | 824 | com = sqrt(P->one_es); |
200 | 824 | if (fabs(P->phi0) > EPS) { |
201 | 174 | sinph0 = sin(P->phi0); |
202 | 174 | cosph0 = cos(P->phi0); |
203 | 174 | con = 1. - P->es * sinph0 * sinph0; |
204 | 174 | Q->B = cosph0 * cosph0; |
205 | 174 | Q->B = sqrt(1. + P->es * Q->B * Q->B / P->one_es); |
206 | 174 | Q->A = Q->B * P->k0 * com / con; |
207 | 174 | D = Q->B * com / (cosph0 * sqrt(con)); |
208 | 174 | if ((F = D * D - 1.) <= 0.) |
209 | 0 | F = 0.; |
210 | 174 | else { |
211 | 174 | F = sqrt(F); |
212 | 174 | if (P->phi0 < 0.) |
213 | 87 | F = -F; |
214 | 174 | } |
215 | 174 | Q->E = F += D; |
216 | 174 | Q->E *= pow(pj_tsfn(P->phi0, sinph0, P->e), Q->B); |
217 | 650 | } else { |
218 | 650 | Q->B = 1. / com; |
219 | 650 | Q->A = P->k0; |
220 | 650 | Q->E = D = F = 1.; |
221 | 650 | } |
222 | 824 | if (alp || gam) { |
223 | 437 | if (alp) { |
224 | 384 | gamma0 = aasin(P->ctx, sin(alpha_c) / D); |
225 | 384 | if (!gam) |
226 | 283 | gamma = alpha_c; |
227 | 384 | } else { |
228 | 53 | gamma0 = gamma; |
229 | 53 | alpha_c = aasin(P->ctx, D * sin(gamma0)); |
230 | 53 | if (proj_errno(P) != 0) { |
231 | | // For a sphere, |gamma| must be <= 90 - |lat_0| |
232 | | // On an ellipsoid, this is very slightly above |
233 | 0 | proj_log_error(P, |
234 | 0 | ("Invalid value for gamma: given lat_0 value, " |
235 | 0 | "|gamma| should be <= " + |
236 | 0 | std::to_string(asin(1. / D) / M_PI * 180) + "°") |
237 | 0 | .c_str()); |
238 | 0 | return pj_default_destructor( |
239 | 0 | P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
240 | 0 | } |
241 | 53 | } |
242 | | |
243 | 437 | if (fabs(fabs(P->phi0) - M_HALFPI) <= TOL) { |
244 | 1 | proj_log_error( |
245 | 1 | P, _("Invalid value for lat_0: |lat_0| should be < 90°")); |
246 | 1 | return pj_default_destructor(P, |
247 | 1 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
248 | 1 | } |
249 | | |
250 | 436 | P->lam0 = lamc - aasin(P->ctx, .5 * (F - 1. / F) * tan(gamma0)) / Q->B; |
251 | 436 | } else { |
252 | 387 | H = pow(pj_tsfn(phi1, sin(phi1), P->e), Q->B); |
253 | 387 | L = pow(pj_tsfn(phi2, sin(phi2), P->e), Q->B); |
254 | 387 | F = Q->E / H; |
255 | 387 | p = (L - H) / (L + H); |
256 | 387 | if (p == 0) { |
257 | | // Not quite, but es is very close to 1... |
258 | 0 | proj_log_error(P, _("Invalid value for eccentricity")); |
259 | 0 | return pj_default_destructor(P, |
260 | 0 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
261 | 0 | } |
262 | 387 | J = Q->E * Q->E; |
263 | 387 | J = (J - L * H) / (J + L * H); |
264 | 387 | if ((con = lam1 - lam2) < -M_PI) |
265 | 20 | lam2 -= M_TWOPI; |
266 | 367 | else if (con > M_PI) |
267 | 9 | lam2 += M_TWOPI; |
268 | 387 | P->lam0 = adjlon(.5 * (lam1 + lam2) - |
269 | 387 | atan(J * tan(.5 * Q->B * (lam1 - lam2)) / p) / Q->B); |
270 | 387 | const double denom = F - 1. / F; |
271 | 387 | if (denom == 0) { |
272 | 0 | proj_log_error(P, _("Invalid value for eccentricity")); |
273 | 0 | return pj_default_destructor(P, |
274 | 0 | PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE); |
275 | 0 | } |
276 | 387 | gamma0 = atan(2. * sin(Q->B * adjlon(lam1 - P->lam0)) / denom); |
277 | 387 | gamma = alpha_c = aasin(P->ctx, D * sin(gamma0)); |
278 | 387 | } |
279 | 823 | Q->singam = sin(gamma0); |
280 | 823 | Q->cosgam = cos(gamma0); |
281 | 823 | Q->sinrot = sin(gamma); |
282 | 823 | Q->cosrot = cos(gamma); |
283 | 823 | Q->BrA = 1. / (Q->ArB = Q->A * (Q->rB = 1. / Q->B)); |
284 | 823 | Q->AB = Q->A * Q->B; |
285 | 823 | if (no_off) |
286 | 74 | Q->u_0 = 0; |
287 | 749 | else { |
288 | 749 | Q->u_0 = fabs(Q->ArB * atan(sqrt(D * D - 1.) / cos(alpha_c))); |
289 | 749 | if (P->phi0 < 0.) |
290 | 87 | Q->u_0 = -Q->u_0; |
291 | 749 | } |
292 | 823 | F = 0.5 * gamma0; |
293 | 823 | Q->v_pole_n = Q->ArB * log(tan(M_FORTPI - F)); |
294 | 823 | Q->v_pole_s = Q->ArB * log(tan(M_FORTPI + F)); |
295 | 823 | P->inv = omerc_e_inverse; |
296 | 823 | P->fwd = omerc_e_forward; |
297 | | |
298 | 823 | return P; |
299 | 824 | } |
300 | | |
301 | | #undef TOL |
302 | | #undef EPS |