Coverage Report

Created: 2025-07-18 07:18

/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
4.82k
#define TOL 1.e-7
44
7.98k
#define EPS 1.e-10
45
46
3.86k
static PJ_XY omerc_e_forward(PJ_LP lp, PJ *P) { /* Ellipsoidal, forward */
47
3.86k
    PJ_XY xy = {0.0, 0.0};
48
3.86k
    struct pj_omerc_data *Q = static_cast<struct pj_omerc_data *>(P->opaque);
49
3.86k
    double u, v;
50
51
3.86k
    if (fabs(fabs(lp.phi) - M_HALFPI) > EPS) {
52
3.86k
        const double W = Q->E / pow(pj_tsfn(lp.phi, sin(lp.phi), P->e), Q->B);
53
3.86k
        const double one_div_W = 1. / W;
54
3.86k
        const double S = .5 * (W - one_div_W);
55
3.86k
        const double T = .5 * (W + one_div_W);
56
3.86k
        const double V = sin(Q->B * lp.lam);
57
3.86k
        const double U = (S * Q->singam - V * Q->cosgam) / T;
58
3.86k
        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
3.86k
        v = 0.5 * Q->ArB * log((1. - U) / (1. + U));
63
3.86k
        const double temp = cos(Q->B * lp.lam);
64
3.86k
        if (fabs(temp) < TOL) {
65
0
            u = Q->A * lp.lam;
66
3.86k
        } else {
67
3.86k
            u = Q->ArB * atan2((S * Q->cosgam + V * Q->singam), temp);
68
3.86k
        }
69
3.86k
    } 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
3.86k
    if (Q->no_rot) {
74
0
        xy.x = u;
75
0
        xy.y = v;
76
3.86k
    } else {
77
3.86k
        u -= Q->u_0;
78
3.86k
        xy.x = v * Q->cosrot + u * Q->sinrot;
79
3.86k
        xy.y = u * Q->cosrot - v * Q->sinrot;
80
3.86k
    }
81
3.86k
    return xy;
82
3.86k
}
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
300
PJ *PJ_PROJECTION(omerc) {
122
300
    double con, com, cosph0, D, F, H, L, sinph0, p, J,
123
300
        gamma = 0, gamma0, lamc = 0, lam1 = 0, lam2 = 0, phi1 = 0, phi2 = 0,
124
300
        alpha_c = 0;
125
300
    int alp, gam, no_off = 0;
126
127
300
    struct pj_omerc_data *Q = static_cast<struct pj_omerc_data *>(
128
300
        calloc(1, sizeof(struct pj_omerc_data)));
129
300
    if (nullptr == Q)
130
0
        return pj_default_destructor(P, PROJ_ERR_OTHER /*ENOMEM*/);
131
300
    P->opaque = Q;
132
133
300
    Q->no_rot = pj_param(P->ctx, P->params, "bno_rot").i;
134
300
    if ((alp = pj_param(P->ctx, P->params, "talpha").i) != 0)
135
91
        alpha_c = pj_param(P->ctx, P->params, "ralpha").f;
136
300
    if ((gam = pj_param(P->ctx, P->params, "tgamma").i) != 0)
137
98
        gamma = pj_param(P->ctx, P->params, "rgamma").f;
138
300
    if (alp || gam) {
139
116
        lamc = pj_param(P->ctx, P->params, "rlonc").f;
140
116
        no_off =
141
            /* For libproj4 compatibility */
142
116
            pj_param(P->ctx, P->params, "tno_off").i
143
            /* for backward compatibility */
144
116
            || pj_param(P->ctx, P->params, "tno_uoff").i;
145
116
        if (no_off) {
146
            /* Mark the parameter as used, so that the pj_get_def() return them
147
             */
148
1
            pj_param(P->ctx, P->params, "sno_uoff");
149
1
            pj_param(P->ctx, P->params, "sno_off");
150
1
        }
151
184
    } else {
152
184
        lam1 = pj_param(P->ctx, P->params, "rlon_1").f;
153
184
        phi1 = pj_param(P->ctx, P->params, "rlat_1").f;
154
184
        lam2 = pj_param(P->ctx, P->params, "rlon_2").f;
155
184
        phi2 = pj_param(P->ctx, P->params, "rlat_2").f;
156
184
        con = fabs(phi1);
157
158
184
        if (fabs(phi1) > M_HALFPI - TOL) {
159
6
            proj_log_error(
160
6
                P, _("Invalid value for lat_1: |lat_1| should be < 90°"));
161
6
            return pj_default_destructor(P,
162
6
                                         PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
163
6
        }
164
178
        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
175
        if (fabs(phi1 - phi2) <= TOL) {
172
20
            proj_log_error(P,
173
20
                           _("Invalid value for lat_1/lat_2: lat_1 should be "
174
20
                             "different from lat_2"));
175
20
            return pj_default_destructor(P,
176
20
                                         PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
177
20
        }
178
179
155
        if (con <= TOL) {
180
6
            proj_log_error(
181
6
                P,
182
6
                _("Invalid value for lat_1: lat_1 should be different from 0"));
183
6
            return pj_default_destructor(P,
184
6
                                         PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
185
6
        }
186
187
149
        if (fabs(fabs(P->phi0) - M_HALFPI) <= TOL) {
188
7
            proj_log_error(
189
7
                P, _("Invalid value for lat_0: |lat_0| should be < 90°"));
190
7
            return pj_default_destructor(P,
191
7
                                         PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
192
7
        }
193
149
    }
194
195
258
    if (pj_param(P->ctx, P->params, "rlon_0").i) {
196
11
        proj_log_trace(P, _("lon_0 is ignored."));
197
11
    }
198
199
258
    com = sqrt(P->one_es);
200
258
    if (fabs(P->phi0) > EPS) {
201
18
        sinph0 = sin(P->phi0);
202
18
        cosph0 = cos(P->phi0);
203
18
        con = 1. - P->es * sinph0 * sinph0;
204
18
        Q->B = cosph0 * cosph0;
205
18
        Q->B = sqrt(1. + P->es * Q->B * Q->B / P->one_es);
206
18
        Q->A = Q->B * P->k0 * com / con;
207
18
        D = Q->B * com / (cosph0 * sqrt(con));
208
18
        if ((F = D * D - 1.) <= 0.)
209
0
            F = 0.;
210
18
        else {
211
18
            F = sqrt(F);
212
18
            if (P->phi0 < 0.)
213
0
                F = -F;
214
18
        }
215
18
        Q->E = F += D;
216
18
        Q->E *= pow(pj_tsfn(P->phi0, sinph0, P->e), Q->B);
217
240
    } else {
218
240
        Q->B = 1. / com;
219
240
        Q->A = P->k0;
220
240
        Q->E = D = F = 1.;
221
240
    }
222
258
    if (alp || gam) {
223
116
        if (alp) {
224
91
            gamma0 = aasin(P->ctx, sin(alpha_c) / D);
225
91
            if (!gam)
226
18
                gamma = alpha_c;
227
91
        } else {
228
25
            gamma0 = gamma;
229
25
            alpha_c = aasin(P->ctx, D * sin(gamma0));
230
25
            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
25
        }
242
243
116
        if (fabs(fabs(P->phi0) - M_HALFPI) <= TOL) {
244
0
            proj_log_error(
245
0
                P, _("Invalid value for lat_0: |lat_0| should be < 90°"));
246
0
            return pj_default_destructor(P,
247
0
                                         PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
248
0
        }
249
250
116
        P->lam0 = lamc - aasin(P->ctx, .5 * (F - 1. / F) * tan(gamma0)) / Q->B;
251
142
    } else {
252
142
        H = pow(pj_tsfn(phi1, sin(phi1), P->e), Q->B);
253
142
        L = pow(pj_tsfn(phi2, sin(phi2), P->e), Q->B);
254
142
        F = Q->E / H;
255
142
        p = (L - H) / (L + H);
256
142
        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
142
        J = Q->E * Q->E;
263
142
        J = (J - L * H) / (J + L * H);
264
142
        if ((con = lam1 - lam2) < -M_PI)
265
7
            lam2 -= M_TWOPI;
266
135
        else if (con > M_PI)
267
23
            lam2 += M_TWOPI;
268
142
        P->lam0 = adjlon(.5 * (lam1 + lam2) -
269
142
                         atan(J * tan(.5 * Q->B * (lam1 - lam2)) / p) / Q->B);
270
142
        const double denom = F - 1. / F;
271
142
        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
142
        gamma0 = atan(2. * sin(Q->B * adjlon(lam1 - P->lam0)) / denom);
277
142
        gamma = alpha_c = aasin(P->ctx, D * sin(gamma0));
278
142
    }
279
258
    Q->singam = sin(gamma0);
280
258
    Q->cosgam = cos(gamma0);
281
258
    Q->sinrot = sin(gamma);
282
258
    Q->cosrot = cos(gamma);
283
258
    Q->BrA = 1. / (Q->ArB = Q->A * (Q->rB = 1. / Q->B));
284
258
    Q->AB = Q->A * Q->B;
285
258
    if (no_off)
286
1
        Q->u_0 = 0;
287
257
    else {
288
257
        Q->u_0 = fabs(Q->ArB * atan(sqrt(D * D - 1.) / cos(alpha_c)));
289
257
        if (P->phi0 < 0.)
290
0
            Q->u_0 = -Q->u_0;
291
257
    }
292
258
    F = 0.5 * gamma0;
293
258
    Q->v_pole_n = Q->ArB * log(tan(M_FORTPI - F));
294
258
    Q->v_pole_s = Q->ArB * log(tan(M_FORTPI + F));
295
258
    P->inv = omerc_e_inverse;
296
258
    P->fwd = omerc_e_forward;
297
298
258
    return P;
299
258
}
300
301
#undef TOL
302
#undef EPS