Coverage Report

Created: 2025-12-31 06:48

next uncovered line (L), next uncovered region (R), next uncovered branch (B)
/src/proj/src/projections/stere.cpp
Line
Count
Source
1
2
#include "proj.h"
3
#include "proj_internal.h"
4
#include <errno.h>
5
#include <math.h>
6
7
PROJ_HEAD(stere, "Stereographic") "\n\tAzi, Sph&Ell\n\tlat_ts=";
8
PROJ_HEAD(ups, "Universal Polar Stereographic") "\n\tAzi, Ell\n\tsouth";
9
10
namespace { // anonymous namespace
11
enum Mode { S_POLE = 0, N_POLE = 1, OBLIQ = 2, EQUIT = 3 };
12
} // anonymous namespace
13
14
namespace { // anonymous namespace
15
struct pj_stere {
16
    double phits;
17
    double sinX1;
18
    double cosX1;
19
    double akm1;
20
    enum Mode mode;
21
};
22
} // anonymous namespace
23
24
3
#define sinph0 static_cast<struct pj_stere *>(P->opaque)->sinX1
25
3
#define cosph0 static_cast<struct pj_stere *>(P->opaque)->cosX1
26
14.7k
#define EPS10 1.e-10
27
0
#define TOL 1.e-8
28
0
#define NITER 8
29
0
#define CONV 1.e-10
30
31
373
static double ssfn_(double phit, double sinphi, double eccen) {
32
373
    sinphi *= eccen;
33
373
    return (tan(.5 * (M_HALFPI + phit)) *
34
373
            pow((1. - sinphi) / (1. + sinphi), .5 * eccen));
35
373
}
36
37
0
static PJ_XY stere_e_forward(PJ_LP lp, PJ *P) { /* Ellipsoidal, forward */
38
0
    PJ_XY xy = {0.0, 0.0};
39
0
    struct pj_stere *Q = static_cast<struct pj_stere *>(P->opaque);
40
0
    double coslam, sinlam, sinX = 0.0, cosX = 0.0, A = 0.0, sinphi;
41
42
0
    coslam = cos(lp.lam);
43
0
    sinlam = sin(lp.lam);
44
0
    sinphi = sin(lp.phi);
45
0
    if (Q->mode == OBLIQ || Q->mode == EQUIT) {
46
0
        const double X = 2. * atan(ssfn_(lp.phi, sinphi, P->e)) - M_HALFPI;
47
0
        sinX = sin(X);
48
0
        cosX = cos(X);
49
0
    }
50
51
0
    switch (Q->mode) {
52
0
    case OBLIQ: {
53
0
        const double denom =
54
0
            Q->cosX1 * (1. + Q->sinX1 * sinX + Q->cosX1 * cosX * coslam);
55
0
        if (denom == 0) {
56
0
            proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN);
57
0
            return proj_coord_error().xy;
58
0
        }
59
0
        A = Q->akm1 / denom;
60
0
        xy.y = A * (Q->cosX1 * sinX - Q->sinX1 * cosX * coslam);
61
0
        xy.x = A * cosX;
62
0
        break;
63
0
    }
64
65
0
    case EQUIT:
66
        /* avoid zero division */
67
0
        if (1. + cosX * coslam == 0.0) {
68
0
            xy.y = HUGE_VAL;
69
0
        } else {
70
0
            A = Q->akm1 / (1. + cosX * coslam);
71
0
            xy.y = A * sinX;
72
0
        }
73
0
        xy.x = A * cosX;
74
0
        break;
75
76
0
    case S_POLE:
77
0
        lp.phi = -lp.phi;
78
0
        coslam = -coslam;
79
0
        sinphi = -sinphi;
80
0
        PROJ_FALLTHROUGH;
81
0
    case N_POLE:
82
0
        if (fabs(lp.phi - M_HALFPI) < 1e-15)
83
0
            xy.x = 0;
84
0
        else
85
0
            xy.x = Q->akm1 * pj_tsfn(lp.phi, sinphi, P->e);
86
0
        xy.y = -xy.x * coslam;
87
0
        break;
88
0
    }
89
90
0
    xy.x = xy.x * sinlam;
91
0
    return xy;
92
0
}
93
94
0
static PJ_XY stere_s_forward(PJ_LP lp, PJ *P) { /* Spheroidal, forward */
95
0
    PJ_XY xy = {0.0, 0.0};
96
0
    struct pj_stere *Q = static_cast<struct pj_stere *>(P->opaque);
97
0
    double sinphi, cosphi, coslam, sinlam;
98
99
0
    sinphi = sin(lp.phi);
100
0
    cosphi = cos(lp.phi);
101
0
    coslam = cos(lp.lam);
102
0
    sinlam = sin(lp.lam);
103
104
0
    switch (Q->mode) {
105
0
    case EQUIT:
106
0
        xy.y = 1. + cosphi * coslam;
107
0
        goto oblcon;
108
0
    case OBLIQ:
109
0
        xy.y = 1. + sinph0 * sinphi + cosph0 * cosphi * coslam;
110
0
    oblcon:
111
0
        if (xy.y <= EPS10) {
112
0
            proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN);
113
0
            return xy;
114
0
        }
115
0
        xy.y = Q->akm1 / xy.y;
116
0
        xy.x = xy.y * cosphi * sinlam;
117
0
        xy.y *= (Q->mode == EQUIT) ? sinphi
118
0
                                   : cosph0 * sinphi - sinph0 * cosphi * coslam;
119
0
        break;
120
0
    case N_POLE:
121
0
        coslam = -coslam;
122
0
        lp.phi = -lp.phi;
123
0
        PROJ_FALLTHROUGH;
124
0
    case S_POLE:
125
0
        if (fabs(lp.phi - M_HALFPI) < TOL) {
126
0
            proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN);
127
0
            return xy;
128
0
        }
129
0
        xy.y = Q->akm1 * tan(M_FORTPI + .5 * lp.phi);
130
0
        xy.x = sinlam * xy.y;
131
0
        xy.y *= coslam;
132
0
        break;
133
0
    }
134
0
    return xy;
135
0
}
136
137
0
static PJ_LP stere_e_inverse(PJ_XY xy, PJ *P) { /* Ellipsoidal, inverse */
138
0
    PJ_LP lp = {0.0, 0.0};
139
0
    struct pj_stere *Q = static_cast<struct pj_stere *>(P->opaque);
140
0
    double cosphi, sinphi, tp = 0.0, phi_l = 0.0, rho, halfe = 0.0,
141
0
                           halfpi = 0.0;
142
143
0
    rho = hypot(xy.x, xy.y);
144
145
0
    switch (Q->mode) {
146
0
    case OBLIQ:
147
0
    case EQUIT:
148
0
        tp = 2. * atan2(rho * Q->cosX1, Q->akm1);
149
0
        cosphi = cos(tp);
150
0
        sinphi = sin(tp);
151
0
        if (rho == 0.0)
152
0
            phi_l = asin(cosphi * Q->sinX1);
153
0
        else
154
0
            phi_l = asin(cosphi * Q->sinX1 + (xy.y * sinphi * Q->cosX1 / rho));
155
156
0
        tp = tan(.5 * (M_HALFPI + phi_l));
157
0
        xy.x *= sinphi;
158
0
        xy.y = rho * Q->cosX1 * cosphi - xy.y * Q->sinX1 * sinphi;
159
0
        halfpi = M_HALFPI;
160
0
        halfe = .5 * P->e;
161
0
        break;
162
0
    case N_POLE:
163
0
        xy.y = -xy.y;
164
0
        PROJ_FALLTHROUGH;
165
0
    case S_POLE:
166
0
        tp = -rho / Q->akm1;
167
0
        phi_l = M_HALFPI - 2. * atan(tp);
168
0
        halfpi = -M_HALFPI;
169
0
        halfe = -.5 * P->e;
170
0
        break;
171
0
    }
172
173
0
    for (int i = NITER; i > 0; --i) {
174
0
        sinphi = P->e * sin(phi_l);
175
0
        lp.phi =
176
0
            2. * atan(tp * pow((1. + sinphi) / (1. - sinphi), halfe)) - halfpi;
177
0
        if (fabs(phi_l - lp.phi) < CONV) {
178
0
            if (Q->mode == S_POLE)
179
0
                lp.phi = -lp.phi;
180
0
            lp.lam = (xy.x == 0. && xy.y == 0.) ? 0. : atan2(xy.x, xy.y);
181
0
            return lp;
182
0
        }
183
0
        phi_l = lp.phi;
184
0
    }
185
186
0
    proj_errno_set(P, PROJ_ERR_COORD_TRANSFM_OUTSIDE_PROJECTION_DOMAIN);
187
0
    return lp;
188
0
}
189
190
0
static PJ_LP stere_s_inverse(PJ_XY xy, PJ *P) { /* Spheroidal, inverse */
191
0
    PJ_LP lp = {0.0, 0.0};
192
0
    struct pj_stere *Q = static_cast<struct pj_stere *>(P->opaque);
193
0
    double c, sinc, cosc;
194
195
0
    const double rh = hypot(xy.x, xy.y);
196
0
    c = 2. * atan(rh / Q->akm1);
197
0
    sinc = sin(c);
198
0
    cosc = cos(c);
199
0
    lp.lam = 0.;
200
201
0
    switch (Q->mode) {
202
0
    case EQUIT:
203
0
        if (fabs(rh) <= EPS10)
204
0
            lp.phi = 0.;
205
0
        else
206
0
            lp.phi = asin(xy.y * sinc / rh);
207
0
        if (cosc != 0. || xy.x != 0.)
208
0
            lp.lam = atan2(xy.x * sinc, cosc * rh);
209
0
        break;
210
0
    case OBLIQ:
211
0
        if (fabs(rh) <= EPS10)
212
0
            lp.phi = P->phi0;
213
0
        else
214
0
            lp.phi = asin(cosc * sinph0 + xy.y * sinc * cosph0 / rh);
215
0
        c = cosc - sinph0 * sin(lp.phi);
216
0
        if (c != 0. || xy.x != 0.)
217
0
            lp.lam = atan2(xy.x * sinc * cosph0, c * rh);
218
0
        break;
219
0
    case N_POLE:
220
0
        xy.y = -xy.y;
221
0
        PROJ_FALLTHROUGH;
222
0
    case S_POLE:
223
0
        if (fabs(rh) <= EPS10)
224
0
            lp.phi = P->phi0;
225
0
        else
226
0
            lp.phi = asin(Q->mode == S_POLE ? -cosc : cosc);
227
0
        lp.lam = (xy.x == 0. && xy.y == 0.) ? 0. : atan2(xy.x, xy.y);
228
0
        break;
229
0
    }
230
0
    return lp;
231
0
}
232
233
7.38k
static PJ *stere_setup(PJ *P) { /* general initialization */
234
7.38k
    double t;
235
7.38k
    struct pj_stere *Q = static_cast<struct pj_stere *>(P->opaque);
236
237
7.38k
    if (fabs((t = fabs(P->phi0)) - M_HALFPI) < EPS10)
238
6.75k
        Q->mode = P->phi0 < 0. ? S_POLE : N_POLE;
239
630
    else
240
630
        Q->mode = t > EPS10 ? OBLIQ : EQUIT;
241
7.38k
    Q->phits = fabs(Q->phits);
242
243
7.38k
    if (P->es != 0.0) {
244
7.06k
        double X;
245
246
7.06k
        switch (Q->mode) {
247
6.59k
        case N_POLE:
248
6.69k
        case S_POLE:
249
6.69k
            if (fabs(Q->phits - M_HALFPI) < EPS10)
250
6.67k
                Q->akm1 =
251
6.67k
                    2. * P->k0 /
252
6.67k
                    sqrt(pow(1 + P->e, 1 + P->e) * pow(1 - P->e, 1 - P->e));
253
23
            else {
254
23
                t = sin(Q->phits);
255
23
                Q->akm1 = cos(Q->phits) / pj_tsfn(Q->phits, t, P->e);
256
23
                t *= P->e;
257
23
                Q->akm1 /= sqrt(1. - t * t);
258
23
            }
259
6.69k
            break;
260
370
        case EQUIT:
261
373
        case OBLIQ:
262
373
            t = sin(P->phi0);
263
373
            X = 2. * atan(ssfn_(P->phi0, t, P->e)) - M_HALFPI;
264
373
            t *= P->e;
265
373
            Q->akm1 = 2. * P->k0 * cos(P->phi0) / sqrt(1. - t * t);
266
373
            Q->sinX1 = sin(X);
267
373
            Q->cosX1 = cos(X);
268
373
            break;
269
7.06k
        }
270
7.06k
        P->inv = stere_e_inverse;
271
7.06k
        P->fwd = stere_e_forward;
272
7.06k
    } else {
273
314
        switch (Q->mode) {
274
3
        case OBLIQ:
275
3
            sinph0 = sin(P->phi0);
276
3
            cosph0 = cos(P->phi0);
277
3
            PROJ_FALLTHROUGH;
278
257
        case EQUIT:
279
257
            Q->akm1 = 2. * P->k0;
280
257
            break;
281
8
        case S_POLE:
282
57
        case N_POLE:
283
57
            Q->akm1 = fabs(Q->phits - M_HALFPI) >= EPS10
284
57
                          ? cos(Q->phits) / tan(M_FORTPI - .5 * Q->phits)
285
57
                          : 2. * P->k0;
286
57
            break;
287
314
        }
288
289
314
        P->inv = stere_s_inverse;
290
314
        P->fwd = stere_s_forward;
291
314
    }
292
7.38k
    return P;
293
7.38k
}
294
295
714
PJ *PJ_PROJECTION(stere) {
296
714
    struct pj_stere *Q =
297
714
        static_cast<struct pj_stere *>(calloc(1, sizeof(struct pj_stere)));
298
714
    if (nullptr == Q)
299
0
        return pj_default_destructor(P, PROJ_ERR_OTHER /*ENOMEM*/);
300
714
    P->opaque = Q;
301
302
714
    Q->phits = pj_param(P->ctx, P->params, "tlat_ts").i
303
714
                   ? pj_param(P->ctx, P->params, "rlat_ts").f
304
714
                   : M_HALFPI;
305
306
714
    return stere_setup(P);
307
714
}
308
309
6.67k
PJ *PJ_PROJECTION(ups) {
310
6.67k
    struct pj_stere *Q =
311
6.67k
        static_cast<struct pj_stere *>(calloc(1, sizeof(struct pj_stere)));
312
6.67k
    if (nullptr == Q)
313
0
        return pj_default_destructor(P, PROJ_ERR_OTHER /*ENOMEM*/);
314
6.67k
    P->opaque = Q;
315
316
    /* International Ellipsoid */
317
6.67k
    P->phi0 = pj_param(P->ctx, P->params, "bsouth").i ? -M_HALFPI : M_HALFPI;
318
6.67k
    if (P->es == 0.0) {
319
5
        proj_log_error(
320
5
            P,
321
5
            _("Invalid value for es: only ellipsoidal formulation supported"));
322
5
        return pj_default_destructor(P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
323
5
    }
324
6.66k
    P->k0 = .994;
325
6.66k
    P->x0 = 2000000.;
326
6.66k
    P->y0 = 2000000.;
327
6.66k
    Q->phits = M_HALFPI;
328
6.66k
    P->lam0 = 0.;
329
330
6.66k
    return stere_setup(P);
331
6.67k
}
332
333
#undef sinph0
334
#undef cosph0
335
#undef EPS10
336
#undef TOL
337
#undef NITER
338
#undef CONV