aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--docs/source/operations/projections/lcc.rst51
-rw-r--r--src/PJ_lcc.c33
-rw-r--r--test/gie/builtins.gie91
-rw-r--r--test/gie/more_builtins.gie17
4 files changed, 172 insertions, 20 deletions
diff --git a/docs/source/operations/projections/lcc.rst b/docs/source/operations/projections/lcc.rst
index 6db2a73c..3b2032d6 100644
--- a/docs/source/operations/projections/lcc.rst
+++ b/docs/source/operations/projections/lcc.rst
@@ -4,6 +4,39 @@
Lambert Conformal Conic
********************************************************************************
+A Lambert Conformal Conic projection (LCC) is a conic map projection
+used for aeronautical charts, portions of the State Plane Coordinate
+System, and many national and regional mapping systems. It is one of
+seven projections introduced by Johann Heinrich Lambert in 1772.
+
+It has several different forms: with one and two standard parallels
+(referred to as 1SP and 2SP in EPSG guidance notes). Additionally we
+provide "2SP Michigan" form which is very similar to normal 2SP, but
+with a scaling factor on the ellipsoid (given as `k_0` parameter).
+It is implemented as per EPSG Guidance Note 7-2 (version 54, August
+2018, page 25). It is used in a few systems in the EPSG database which
+justifies adding this otherwise non-standard projection.
+
++---------------------+----------------------------------------------------------+
+| **Classification** | Conformal conic |
++---------------------+----------------------------------------------------------+
+| **Available forms** | Forward and inverse, spherical and elliptical projection.|
+| | One or two standard parallels (1SP and 2SP). |
+| | "LCC 2SP Michigan" form can be used by setting |
+| | `+k_0` parameter to specify elliposid scale. |
++---------------------+----------------------------------------------------------+
+| **Defined area** | Best for regions predominantly east–west in extent and |
+| | located in the middle north or south latitudes. |
++---------------------+----------------------------------------------------------+
+| **Alias** | lcc |
++---------------------+----------------------------------------------------------+
+| **Domain** | 2D |
++---------------------+----------------------------------------------------------+
+| **Input type** | Geodetic coordinates |
++---------------------+----------------------------------------------------------+
+| **Output type** | Projected coordinates |
++---------------------+----------------------------------------------------------+
+
.. figure:: ./images/lcc.png
:width: 500 px
:align: center
@@ -31,3 +64,21 @@ Parameters
.. include:: ../options/x_0.rst
.. include:: ../options/y_0.rst
+
+.. option:: +k_0=<value>
+
+ This parameter can represent two different values depending on the
+ form of the projection. In LCC 1SP it determines the scale factor
+ at natural origin. In LCC 2SP Michigan it determines the ellipsoid
+ scale factor.
+
+ *Defaults to 1.0.*
+
+Further reading
+###############
+
+#. `Wikipedia <https://en.wikipedia.org/wiki/Lambert_conformal_conic_projection>`_
+#. `Wolfram Mathworld <http://mathworld.wolfram.com/LambertConformalConicProjection.html>`_
+#. `John P. Snyder "Map projections: A working manual" (pp. 104-110) <https://pubs.er.usgs.gov/publication/pp1395>`_
+#. `ArcGIS documentation on "Lambert Conformal Conic" <http://desktop.arcgis.com/en/arcmap/10.3/guide-books/map-projections/lambert-conformal-conic.htm>`_
+#. `EPSG Guidance Note 7-2 (version 54, August 2018, page 25) <http://www.epsg.org/Guidancenotes.aspx>`_
diff --git a/src/PJ_lcc.c b/src/PJ_lcc.c
index 34ed99cb..96aa3f2a 100644
--- a/src/PJ_lcc.c
+++ b/src/PJ_lcc.c
@@ -5,9 +5,9 @@
#include "proj_math.h"
PROJ_HEAD(lcc, "Lambert Conformal Conic")
- "\n\tConic, Sph&Ell\n\tlat_1= and lat_2= or lat_0";
+ "\n\tConic, Sph&Ell\n\tlat_1= and lat_2= or lat_0, k_0=";
-# define EPS10 1.e-10
+#define EPS10 1.e-10
struct pj_opaque {
double phi1;
@@ -15,12 +15,11 @@ struct pj_opaque {
double n;
double rho0;
double c;
- int ellips;
};
static XY e_forward (LP lp, PJ *P) { /* Ellipsoidal, forward */
- XY xy = {0.0,0.0};
+ XY xy = {0., 0.};
struct pj_opaque *Q = P->opaque;
double rho;
@@ -31,18 +30,19 @@ static XY e_forward (LP lp, PJ *P) { /* Ellipsoidal, forward */
}
rho = 0.;
} else {
- rho = Q->c * (Q->ellips ? pow(pj_tsfn(lp.phi, sin(lp.phi),
- P->e), Q->n) : pow(tan(M_FORTPI + .5 * lp.phi), -Q->n));
+ rho = Q->c * (P->es != 0. ?
+ pow(pj_tsfn(lp.phi, sin(lp.phi), P->e), Q->n) :
+ pow(tan(M_FORTPI + .5 * lp.phi), -Q->n));
}
lp.lam *= Q->n;
- xy.x = P->k0 * (rho * sin( lp.lam) );
- xy.y = P->k0 * (Q->rho0 - rho * cos(lp.lam) );
+ xy.x = P->k0 * (rho * sin(lp.lam));
+ xy.y = P->k0 * (Q->rho0 - rho * cos(lp.lam));
return xy;
}
static LP e_inverse (XY xy, PJ *P) { /* Ellipsoidal, inverse */
- LP lp = {0.0,0.0};
+ LP lp = {0., 0.};
struct pj_opaque *Q = P->opaque;
double rho;
@@ -51,13 +51,13 @@ static LP e_inverse (XY xy, PJ *P) { /* Ellipsoidal, inverse */
xy.y = Q->rho0 - xy.y;
rho = hypot(xy.x, xy.y);
- if (rho != 0.0) {
+ if (rho != 0.) {
if (Q->n < 0.) {
rho = -rho;
xy.x = -xy.x;
xy.y = -xy.y;
}
- if (Q->ellips) {
+ if (P->es != 0.) {
lp.phi = pj_phi2(P->ctx, pow(rho / Q->c, 1./Q->n), P->e);
if (lp.phi == HUGE_VAL) {
proj_errno_set(P, PJD_ERR_TOLERANCE_CONDITION);
@@ -78,13 +78,12 @@ static LP e_inverse (XY xy, PJ *P) { /* Ellipsoidal, inverse */
PJ *PROJECTION(lcc) {
double cosphi, sinphi;
int secant;
- struct pj_opaque *Q = pj_calloc (1, sizeof (struct pj_opaque));
+ struct pj_opaque *Q = pj_calloc(1, sizeof (struct pj_opaque));
- if (0==Q)
- return pj_default_destructor (P, ENOMEM);
+ if (0 == Q)
+ return pj_default_destructor(P, ENOMEM);
P->opaque = Q;
-
Q->phi1 = pj_param(P->ctx, P->params, "rlat_1").f;
if (pj_param(P->ctx, P->params, "tlat_2").i)
Q->phi2 = pj_param(P->ctx, P->params, "rlat_2").f;
@@ -99,10 +98,9 @@ PJ *PROJECTION(lcc) {
Q->n = sinphi = sin(Q->phi1);
cosphi = cos(Q->phi1);
secant = fabs(Q->phi1 - Q->phi2) >= EPS10;
- if( (Q->ellips = (P->es != 0.)) ) {
+ if (P->es != 0.) {
double ml1, m1;
- P->e = sqrt(P->es);
m1 = pj_msfn(sinphi, cosphi, P->es);
ml1 = pj_tsfn(Q->phi1, sinphi, P->e);
if (secant) { /* secant cone */
@@ -128,4 +126,3 @@ PJ *PROJECTION(lcc) {
return P;
}
-
diff --git a/test/gie/builtins.gie b/test/gie/builtins.gie
index 1ee968f4..955b5b51 100644
--- a/test/gie/builtins.gie
+++ b/test/gie/builtins.gie
@@ -2505,11 +2505,10 @@ expect -217928.275907355 112144.329220142
accept -2 -1
expect -217928.275907355 -112144.329220142
-
===============================================================================
Lambert Conformal Conic
Conic, Sph&Ell
- lat_1= and lat_2= or lat_0
+ lat_1= and lat_2= or lat_0, k_0=
===============================================================================
-------------------------------------------------------------------------------
@@ -2535,6 +2534,94 @@ expect -0.001796359 0.000904232
accept -200 -100
expect -0.001796358 -0.000904233
+-------------------------------------------------------------------------------
+Tests with +k_0 (Lambert Conformal Conic 2SP Michigan).
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=GRS80 +lat_1=0.5 +lat_2=2 +k_0=1.0000382
+-------------------------------------------------------------------------------
+tolerance 0.1 mm
+accept 2 1
+expect 222596.942614366 110664.761103214
+accept 2 -1
+expect 222765.389013083 -110537.020013748
+accept -2 1
+expect -222596.942614366 110664.761103214
+accept -2 -1
+expect -222765.389013083 -110537.020013748
+
+direction inverse
+accept 200 100
+expect 0.001796291 0.000904198
+accept 200 -100
+expect 0.001796290 -0.000904199
+accept -200 100
+expect -0.001796291 0.000904198
+accept -200 -100
+expect -0.001796290 -0.000904199
+
+
+-------------------------------------------------------------------------------
+Test various corner cases, spherical projection and one standard
+parallel to improve test coverage.
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=GRS80 +lat_1=0.5 +lat_2=2
+-------------------------------------------------------------------------------
+tolerance 0.1 mm
+
+accept 0 90
+expect 0 292411117.537843227
+accept 0 0
+expect 0 0
+
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=sphere +lat_1=30 +lat_2=40
+-------------------------------------------------------------------------------
+tolerance 0.1 mm
+
+accept 1 2
+expect 129391.909521100 262101.674176860
+accept 0 0
+expect 0 0
+
+direction inverse
+
+accept 129391.909521100 262101.674176860
+expect 1 2
+accept 0 0
+expect 0 0
+
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=GRS80 +lat_1=30
+-------------------------------------------------------------------------------
+tolerance 0.1 mm
+
+accept 1 2
+expect 131833.493971117 265456.213515346
+accept 1 -2
+expect 137536.205750651 -269686.591917190
+accept -1 2
+expect -131833.493971117 265456.213515346
+accept -1 -2
+expect -137536.205750651 -269686.591917190
+
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=sphere +lat_1=30
+-------------------------------------------------------------------------------
+tolerance 0.1 mm
+
+accept 1 2
+expect 131824.206082557 267239.875053699
+accept -1 -2
+expect -137565.475967350 -271546.945608449
+accept 1 -2
+expect 137565.475967350 -271546.945608449
+accept -1 2
+expect -131824.206082557 267239.875053699
+
+direction inverse
+
+accept 131824.206082557 267239.875053699
+expect 1 2
===============================================================================
Lambert Conformal Conic Alternative
diff --git a/test/gie/more_builtins.gie b/test/gie/more_builtins.gie
index 10feea70..960e5814 100644
--- a/test/gie/more_builtins.gie
+++ b/test/gie/more_builtins.gie
@@ -243,6 +243,23 @@ expect failure errno no_args
+-------------------------------------------------------------------------------
+Tests for LCC 2SP Michigan (from PJ_lcc.c)
+-------------------------------------------------------------------------------
+This test is taken from EPSG guidance note 7-2 (version 54, August 2018,
+page 25)
+-------------------------------------------------------------------------------
+operation +proj=lcc +ellps=clrk66 +lat_1=44d11'N +lat_2=45d42'N +x_0=609601.2192 +lon_0=84d20'W +lat_0=43d19'N +k_0=1.0000382 +units=us-ft
+-------------------------------------------------------------------------------
+tolerance 5 mm
+accept 83d10'W 43d45'N
+expect 2308335.75 160210.48
+
+direction inverse
+accept 2308335.75 160210.48
+expect 83d10'W 43d45'N
+-------------------------------------------------------------------------------
+
-------------------------------------------------------------------------------
A number of tests from PJ_helmert.c