-
Notifications
You must be signed in to change notification settings - Fork 362
/
crs.py
3234 lines (2686 loc) · 113 KB
/
crs.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# Copyright Crown and Cartopy Contributors
#
# This file is part of Cartopy and is released under the BSD 3-clause license.
# See LICENSE in the root of the repository for full licensing details.
"""
The crs module defines Coordinate Reference Systems and the transformations
between them.
"""
from abc import ABCMeta
from collections import OrderedDict
from functools import lru_cache
import io
import json
import math
import warnings
import numpy as np
from pyproj import Transformer
from pyproj.exceptions import ProjError
import shapely.geometry as sgeom
from shapely.prepared import prep
import cartopy.trace
try:
# https://github.com/pyproj4/pyproj/pull/912
from pyproj.crs import CustomConstructorCRS as _CRS
except ImportError:
from pyproj import CRS as _CRS
__document_these__ = ['CRS', 'Geocentric', 'Geodetic', 'Globe']
WGS84_SEMIMAJOR_AXIS = 6378137.0
WGS84_SEMIMINOR_AXIS = 6356752.3142
# Cache the transformer creation method
@lru_cache
def _get_transformer_from_crs(src_crs, tgt_crs):
return Transformer.from_crs(src_crs, tgt_crs, always_xy=True)
def _safe_pj_transform(src_crs, tgt_crs, x, y, z=None, trap=True):
transformer = _get_transformer_from_crs(src_crs, tgt_crs)
# if a projection is essentially 2d there
# should be no harm in setting its z to 0
if z is None:
z = np.zeros_like(x)
with warnings.catch_warnings():
# pyproj implicitly converts size-1 arrays to scalars, which is
# deprecated in numpy 1.25, but *also* handles the future error
# see https://github.com/numpy/numpy/pull/10615
# and https://github.com/SciTools/cartopy/pull/2194
warnings.filterwarnings(
"ignore",
message="Conversion of an array with ndim > 0"
)
return transformer.transform(x, y, z, errcheck=trap)
class Globe:
"""
Define an ellipsoid and, optionally, how to relate it to the real world.
"""
def __init__(self, datum=None, ellipse='WGS84',
semimajor_axis=None, semiminor_axis=None,
flattening=None, inverse_flattening=None,
towgs84=None, nadgrids=None):
"""
Parameters
----------
datum
Proj "datum" definition. Defaults to None.
ellipse
Proj "ellps" definition. Defaults to 'WGS84'.
semimajor_axis
Semimajor axis of the spheroid / ellipsoid. Defaults to None.
semiminor_axis
Semiminor axis of the ellipsoid. Defaults to None.
flattening
Flattening of the ellipsoid. Defaults to None.
inverse_flattening
Inverse flattening of the ellipsoid. Defaults to None.
towgs84
Passed through to the Proj definition. Defaults to None.
nadgrids
Passed through to the Proj definition. Defaults to None.
"""
self.datum = datum
self.ellipse = ellipse
self.semimajor_axis = semimajor_axis
self.semiminor_axis = semiminor_axis
self.flattening = flattening
self.inverse_flattening = inverse_flattening
self.towgs84 = towgs84
self.nadgrids = nadgrids
def to_proj4_params(self):
"""
Create an OrderedDict of key value pairs which represents this globe
in terms of proj params.
"""
proj4_params = (
['datum', self.datum],
['ellps', self.ellipse],
['a', self.semimajor_axis],
['b', self.semiminor_axis],
['f', self.flattening],
['rf', self.inverse_flattening],
['towgs84', self.towgs84],
['nadgrids', self.nadgrids]
)
return OrderedDict((k, v) for k, v in proj4_params if v is not None)
class CRS(_CRS):
"""
Define a Coordinate Reference System using proj. The :class:`cartopy.crs.CRS`
class is the very core of cartopy, all coordinate reference systems in cartopy
have :class:`~cartopy.crs.CRS` as a parent class.
"""
#: Whether this projection can handle ellipses.
_handles_ellipses = True
def __init__(self, proj4_params, globe=None):
"""
Parameters
----------
proj4_params: iterable of key-value pairs
The proj4 parameters required to define the
desired CRS. The parameters should not describe
the desired elliptic model, instead create an
appropriate Globe instance. The ``proj4_params``
parameters will override any parameters that the
Globe defines.
globe: :class:`~cartopy.crs.Globe` instance, optional
If omitted, the default Globe instance will be created.
See :class:`~cartopy.crs.Globe` for details.
"""
self.input = (proj4_params, globe)
# for compatibility with pyproj.CRS and rasterio.crs.CRS
try:
proj4_params = proj4_params.to_wkt()
except AttributeError:
pass
# handle PROJ JSON
if (
isinstance(proj4_params, dict) and
"proj" not in proj4_params and
"init" not in proj4_params
):
proj4_params = json.dumps(proj4_params)
if globe is not None and isinstance(proj4_params, str):
raise ValueError("Cannot have 'globe' with string params.")
if globe is None and not isinstance(proj4_params, str):
if self._handles_ellipses:
globe = Globe()
else:
globe = Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS,
ellipse=None)
if globe is not None and not self._handles_ellipses:
a = globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS
b = globe.semiminor_axis or a
if a != b or globe.ellipse is not None:
warnings.warn(f'The {self.__class__.__name__!r} projection '
'does not handle elliptical globes.')
self.globe = globe
if isinstance(proj4_params, str):
self._proj4_params = {}
self.proj4_init = proj4_params
else:
self._proj4_params = self.globe.to_proj4_params()
self._proj4_params.update(proj4_params)
init_items = []
for k, v in self._proj4_params.items():
if v is not None:
if isinstance(v, float):
init_items.append(f'+{k}={v:.16}')
elif isinstance(v, np.float32):
init_items.append(f'+{k}={v:.8}')
else:
init_items.append(f'+{k}={v}')
else:
init_items.append(f'+{k}')
self.proj4_init = ' '.join(init_items) + ' +no_defs'
super().__init__(self.proj4_init)
def __eq__(self, other):
if isinstance(other, CRS) and self.proj4_init == other.proj4_init:
# Fast path Cartopy's CRS
return True
# For everything else, we let pyproj handle the comparison
return super().__eq__(other)
def __hash__(self):
"""Hash the CRS based on its proj4_init string."""
return hash(self.proj4_init)
def __reduce__(self):
"""
Implement the __reduce__ method used when pickling or performing deepcopy.
"""
if type(self) is CRS:
# State can be reproduced by the proj4_params and globe inputs.
return self.__class__, self.input
else:
# Produces a stateless instance of this class (e.g. an empty tuple).
# The state will then be added via __getstate__ and __setstate__.
# We are forced to this approach because a CRS does not store
# the constructor keyword arguments in its state.
return self.__class__, (), self.__getstate__()
def __getstate__(self):
"""Return the full state of this instance for reconstruction
in ``__setstate__``.
"""
state = self.__dict__.copy()
# remove pyproj specific attrs
state.pop('srs', None)
state.pop('_local', None)
# Remove the proj4 instance and the proj4_init string, which can
# be re-created (in __setstate__) from the other arguments.
state.pop('proj4', None)
state.pop('proj4_init', None)
state['proj4_params'] = self.proj4_params
return state
def __setstate__(self, state):
"""
Take the dictionary created by ``__getstate__`` and passes it
through to this implementation's __init__ method.
"""
# Strip out the key state items for a CRS instance.
CRS_state = {key: state.pop(key) for key in ['proj4_params', 'globe']}
# Put everything else directly into the dict of the instance.
self.__dict__.update(state)
# Call the init of this class to ensure that the projection is
# properly initialised with proj4.
CRS.__init__(self, **CRS_state)
def _as_mpl_transform(self, axes=None):
"""
Cast this CRS instance into a :class:`matplotlib.axes.Axes` using
the Matplotlib ``_as_mpl_transform`` interface.
"""
# lazy import mpl.geoaxes (and therefore matplotlib) as mpl
# is only an optional dependency
import cartopy.mpl.geoaxes as geoaxes
if not isinstance(axes, geoaxes.GeoAxes):
raise ValueError(
'Axes should be an instance of GeoAxes, got %s' % type(axes)
)
return (
geoaxes.InterProjectionTransform(self, axes.projection) +
axes.transData
)
@property
def proj4_params(self):
return dict(self._proj4_params)
def as_geocentric(self):
"""
Return a new Geocentric CRS with the same ellipse/datum as this
CRS.
"""
return CRS(
{
"$schema": (
"https://proj.org/schemas/v0.2/projjson.schema.json"
),
"type": "GeodeticCRS",
"name": "unknown",
"datum": self.datum.to_json_dict(),
"coordinate_system": {
"subtype": "Cartesian",
"axis": [
{
"name": "Geocentric X",
"abbreviation": "X",
"direction": "geocentricX",
"unit": "metre"
},
{
"name": "Geocentric Y",
"abbreviation": "Y",
"direction": "geocentricY",
"unit": "metre"
},
{
"name": "Geocentric Z",
"abbreviation": "Z",
"direction": "geocentricZ",
"unit": "metre"
}
]
}
}
)
def as_geodetic(self):
"""
Return a new Geodetic CRS with the same ellipse/datum as this
CRS.
"""
return CRS(self.geodetic_crs.srs)
def is_geodetic(self):
return self.is_geographic
def transform_point(self, x, y, src_crs, trap=True):
"""
transform_point(x, y, src_crs)
Transform the given float64 coordinate pair, in the given source
coordinate system (``src_crs``), to this coordinate system.
Parameters
----------
x
the x coordinate, in ``src_crs`` coordinates, to transform
y
the y coordinate, in ``src_crs`` coordinates, to transform
src_crs
instance of :class:`CRS` that represents the coordinate
system of ``x`` and ``y``.
trap
Whether proj errors for "latitude or longitude exceeded limits" and
"tolerance condition error" should be trapped.
Returns
-------
(x, y) in this coordinate system
"""
result = self.transform_points(
src_crs, np.array([x]), np.array([y]), trap=trap,
).reshape((1, 3))
return result[0, 0], result[0, 1]
def transform_points(self, src_crs, x, y, z=None, trap=False):
"""
transform_points(src_crs, x, y[, z])
Transform the given coordinates, in the given source
coordinate system (``src_crs``), to this coordinate system.
Parameters
----------
src_crs
instance of :class:`CRS` that represents the
coordinate system of ``x``, ``y`` and ``z``.
x
the x coordinates (array), in ``src_crs`` coordinates,
to transform. May be 1 or 2 dimensional.
y
the y coordinates (array), in ``src_crs`` coordinates,
to transform. Its shape must match that of x.
z: optional
the z coordinates (array), in ``src_crs`` coordinates, to
transform. Defaults to None.
If supplied, its shape must match that of x.
trap
Whether proj errors for "latitude or longitude exceeded limits" and
"tolerance condition error" should be trapped.
Returns
-------
Array of shape ``x.shape + (3, )`` in this coordinate system.
"""
result_shape = tuple(x.shape[i] for i in range(x.ndim)) + (3, )
if z is None:
if x.ndim > 2 or y.ndim > 2:
raise ValueError('x and y arrays must be 1 or 2 dimensional')
elif x.ndim != 1 or y.ndim != 1:
x, y = x.flatten(), y.flatten()
if x.shape[0] != y.shape[0]:
raise ValueError('x and y arrays must have the same length')
else:
if x.ndim > 2 or y.ndim > 2 or z.ndim > 2:
raise ValueError('x, y and z arrays must be 1 or 2 '
'dimensional')
elif x.ndim != 1 or y.ndim != 1 or z.ndim != 1:
x, y, z = x.flatten(), y.flatten(), z.flatten()
if not x.shape[0] == y.shape[0] == z.shape[0]:
raise ValueError('x, y, and z arrays must have the same '
'length')
npts = x.shape[0]
result = np.empty([npts, 3], dtype=np.double)
if npts:
if self == src_crs and (
isinstance(src_crs, _CylindricalProjection) or
self.is_geodetic()):
# convert from [0,360] to [-180,180]
x = np.array(x, copy=True)
to_180 = (x > 180) | (x < -180)
x[to_180] = (((x[to_180] + 180) % 360) - 180)
try:
result[:, 0], result[:, 1], result[:, 2] = \
_safe_pj_transform(src_crs, self, x, y, z, trap=trap)
except ProjError as err:
msg = str(err).lower()
if (
"latitude" in msg or
"longitude" in msg or
"outside of projection domain" in msg or
"tolerance condition error" in msg
):
result[:] = np.nan
else:
raise
if not trap:
result[np.isinf(result)] = np.nan
if len(result_shape) > 2:
return result.reshape(result_shape)
return result
def transform_vectors(self, src_proj, x, y, u, v):
"""
transform_vectors(src_proj, x, y, u, v)
Transform the given vector components, with coordinates in the
given source coordinate system (``src_proj``), to this coordinate
system. The vector components must be given relative to the
source projection's coordinate reference system (grid eastward and
grid northward).
Parameters
----------
src_proj
The :class:`CRS.Projection` that represents the coordinate system
the vectors are defined in.
x
The x coordinates of the vectors in the source projection.
y
The y coordinates of the vectors in the source projection.
u
The grid-eastward components of the vectors.
v
The grid-northward components of the vectors.
Note
----
x, y, u and v may be 1 or 2 dimensional, but must all have matching
shapes.
Returns
-------
ut, vt: The transformed vector components.
Note
----
The algorithm used to transform vectors is an approximation
rather than an exact transform, but the accuracy should be
good enough for visualization purposes.
"""
if not (x.shape == y.shape == u.shape == v.shape):
raise ValueError('x, y, u and v arrays must be the same shape')
if x.ndim not in (1, 2):
raise ValueError('x, y, u and v must be 1 or 2 dimensional')
# Transform the coordinates to the target projection.
proj_xyz = self.transform_points(src_proj, x, y)
target_x, target_y = proj_xyz[..., 0], proj_xyz[..., 1]
# Rotate the input vectors to the projection.
#
# 1: Find the magnitude and direction of the input vectors.
vector_magnitudes = (u**2 + v**2)**0.5
vector_angles = np.arctan2(v, u)
# 2: Find a point in the direction of the original vector that is
# a small distance away from the base point of the vector (near
# the poles the point may have to be in the opposite direction
# to be valid).
factor = 360000.
delta = (src_proj.x_limits[1] - src_proj.x_limits[0]) / factor
x_perturbations = delta * np.cos(vector_angles)
y_perturbations = delta * np.sin(vector_angles)
# 3: Handle points that are invalid. These come from picking a new
# point that is outside the domain of the CRS. The first step is
# to apply the native transform to the input coordinates to make
# sure they are in the appropriate range. Then detect all the
# coordinates where the perturbation takes the point out of the
# valid x-domain and fix them. After that do the same for points
# that are outside the valid y-domain, which may reintroduce some
# points outside of the valid x-domain
proj_xyz = src_proj.transform_points(src_proj, x, y)
source_x, source_y = proj_xyz[..., 0], proj_xyz[..., 1]
# Detect all the coordinates where the perturbation takes the point
# outside of the valid x-domain, and reverse the direction of the
# perturbation to fix this.
eps = 1e-9
invalid_x = np.logical_or(
source_x + x_perturbations < src_proj.x_limits[0] - eps,
source_x + x_perturbations > src_proj.x_limits[1] + eps)
if invalid_x.any():
x_perturbations[invalid_x] *= -1
y_perturbations[invalid_x] *= -1
# Do the same for coordinates where the perturbation takes the point
# outside of the valid y-domain. This may reintroduce some points
# that will be outside the x-domain when the perturbation is
# applied.
invalid_y = np.logical_or(
source_y + y_perturbations < src_proj.y_limits[0] - eps,
source_y + y_perturbations > src_proj.y_limits[1] + eps)
if invalid_y.any():
x_perturbations[invalid_y] *= -1
y_perturbations[invalid_y] *= -1
# Keep track of the points where the perturbation direction was
# reversed.
reversed_vectors = np.logical_xor(invalid_x, invalid_y)
# See if there were any points where we cannot reverse the direction
# of the perturbation to get the perturbed point within the valid
# domain of the projection, and issue a warning if there are.
problem_points = np.logical_or(
source_x + x_perturbations < src_proj.x_limits[0] - eps,
source_x + x_perturbations > src_proj.x_limits[1] + eps)
if problem_points.any():
warnings.warn('Some vectors at source domain corners '
'may not have been transformed correctly')
# 4: Transform this set of points to the projection coordinates and
# find the angle between the base point and the perturbed point
# in the projection coordinates (reversing the direction at any
# points where the original was reversed in step 3).
proj_xyz = self.transform_points(src_proj,
source_x + x_perturbations,
source_y + y_perturbations)
target_x_perturbed = proj_xyz[..., 0]
target_y_perturbed = proj_xyz[..., 1]
projected_angles = np.arctan2(target_y_perturbed - target_y,
target_x_perturbed - target_x)
if reversed_vectors.any():
projected_angles[reversed_vectors] += np.pi
# 5: Form the projected vector components, preserving the magnitude
# of the original vectors.
projected_u = vector_magnitudes * np.cos(projected_angles)
projected_v = vector_magnitudes * np.sin(projected_angles)
return projected_u, projected_v
class Geodetic(CRS):
"""
Define a latitude/longitude coordinate system with spherical topology,
geographical distance and coordinates are measured in degrees.
"""
def __init__(self, globe=None):
"""
Parameters
----------
globe: A :class:`cartopy.crs.Globe`, optional
Defaults to a "WGS84" datum.
"""
proj4_params = [('proj', 'lonlat')]
globe = globe or Globe(datum='WGS84')
super().__init__(proj4_params, globe)
# XXX Implement fwd such as Basemap's Geod.
# Would be used in the tissot example.
# Could come from https://geographiclib.sourceforge.io
class Geocentric(CRS):
"""
Define a Geocentric coordinate system, where x, y, z are Cartesian
coordinates from the center of the Earth.
"""
def __init__(self, globe=None):
"""
Parameters
----------
globe: A :class:`cartopy.crs.Globe`, optional
Defaults to a "WGS84" datum.
"""
proj4_params = [('proj', 'geocent')]
globe = globe or Globe(datum='WGS84')
super().__init__(proj4_params, globe)
class RotatedGeodetic(CRS):
"""
Define a rotated latitude/longitude coordinate system with spherical
topology and geographical distance.
Coordinates are measured in degrees.
The class uses proj to perform an ob_tran operation, using the
pole_longitude to set a lon_0 then performing two rotations based on
pole_latitude and central_rotated_longitude.
This is equivalent to setting the new pole to a location defined by
the pole_latitude and pole_longitude values in the GeogCRS defined by
globe, then rotating this new CRS about it's pole using the
central_rotated_longitude value.
"""
def __init__(self, pole_longitude, pole_latitude,
central_rotated_longitude=0.0, globe=None):
"""
Parameters
----------
pole_longitude
Pole longitude position, in unrotated degrees.
pole_latitude
Pole latitude position, in unrotated degrees.
central_rotated_longitude: optional
Longitude rotation about the new pole, in degrees. Defaults to 0.
globe: optional
A :class:`cartopy.crs.Globe`. Defaults to a "WGS84" datum.
"""
globe = globe or Globe(datum='WGS84')
proj4_params = [('proj', 'ob_tran'), ('o_proj', 'latlon'),
('o_lon_p', central_rotated_longitude),
('o_lat_p', pole_latitude),
('lon_0', 180 + pole_longitude),
('to_meter', math.radians(1) * (
globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS))]
super().__init__(proj4_params, globe=globe)
class Projection(CRS, metaclass=ABCMeta):
"""
Define a projected coordinate system with flat topology and Euclidean
distance.
"""
_method_map = {
'Point': '_project_point',
'LineString': '_project_line_string',
'LinearRing': '_project_linear_ring',
'Polygon': '_project_polygon',
'MultiPoint': '_project_multipoint',
'MultiLineString': '_project_multiline',
'MultiPolygon': '_project_multipolygon',
}
# Whether or not this projection can handle wrapped coordinates
_wrappable = False
def __init__(self, *args, **kwargs):
super().__init__(*args, **kwargs)
self.bounds = None
if self.area_of_use:
# Convert lat/lon bounds to projected bounds.
# Geographic area of the entire dataset referenced to WGS 84
# NB. We can't use a polygon transform at this stage because
# that relies on the existence of the map boundary... the very
# thing we're trying to work out! ;-)
x0 = self.area_of_use.west
x1 = self.area_of_use.east
y0 = self.area_of_use.south
y1 = self.area_of_use.north
lons = np.array([x0, x0, x1, x1])
lats = np.array([y0, y1, y1, y0])
points = self.transform_points(self.as_geodetic(), lons, lats)
x = points[:, 0]
y = points[:, 1]
self.bounds = (x.min(), x.max(), y.min(), y.max())
x0, x1, y0, y1 = self.bounds
self.threshold = min(x1 - x0, y1 - y0) / 100.
@property
def boundary(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return sgeom.LineString([(x0, y0), (x0, y1), (x1, y1), (x1, y0),
(x0, y0)])
@property
def x_limits(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return (x0, x1)
@property
def y_limits(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return (y0, y1)
@property
def threshold(self):
return getattr(self, '_threshold', 0.5)
@threshold.setter
def threshold(self, t):
self._threshold = t
@property
def cw_boundary(self):
try:
boundary = self._cw_boundary
except AttributeError:
boundary = sgeom.LinearRing(self.boundary)
self._cw_boundary = boundary
return boundary
@property
def ccw_boundary(self):
try:
boundary = self._ccw_boundary
except AttributeError:
boundary = sgeom.LinearRing(self.boundary.coords[::-1])
self._ccw_boundary = boundary
return boundary
@property
def domain(self):
try:
domain = self._domain
except AttributeError:
domain = self._domain = sgeom.Polygon(self.boundary)
return domain
def is_geodetic(self):
return False
def _determine_longitude_bounds(self, central_longitude):
# In new proj, using exact limits will wrap-around, so subtract a
# small epsilon:
epsilon = 1e-10
minlon = -180 + central_longitude
maxlon = 180 + central_longitude
if central_longitude > 0:
maxlon -= epsilon
elif central_longitude < 0:
minlon += epsilon
return minlon, maxlon
def _repr_html_(self):
from html import escape
try:
# As matplotlib is not a core cartopy dependency, don't error
# if it's not available.
import matplotlib.pyplot as plt
except ImportError:
# We can't return an SVG of the CRS, so let Jupyter fall back to
# a default repr by returning None.
return None
# Produce a visual repr of the Projection instance.
fig, ax = plt.subplots(figsize=(5, 3),
subplot_kw={'projection': self})
ax.set_global()
ax.coastlines('auto')
ax.gridlines()
buf = io.StringIO()
fig.savefig(buf, format='svg', bbox_inches='tight')
plt.close(fig)
# "Rewind" the buffer to the start and return it as an svg string.
buf.seek(0)
svg = buf.read()
return f'{svg}<pre>{escape(object.__repr__(self))}</pre>'
def _as_mpl_axes(self):
import cartopy.mpl.geoaxes as geoaxes
return geoaxes.GeoAxes, {'projection': self}
def project_geometry(self, geometry, src_crs=None):
"""
Project the given geometry into this projection.
Parameters
----------
geometry
The geometry to (re-)project.
src_crs: optional
The source CRS. Defaults to None.
If src_crs is None, the source CRS is assumed to be a geodetic
version of the target CRS.
Returns
-------
geometry
The projected result (a shapely geometry).
"""
if src_crs is None:
src_crs = self.as_geodetic()
elif not isinstance(src_crs, CRS):
raise TypeError('Source CRS must be an instance of CRS'
' or one of its subclasses, or None.')
geom_type = geometry.geom_type
method_name = self._method_map.get(geom_type)
if not method_name:
raise ValueError(f'Unsupported geometry type {geom_type!r}')
return getattr(self, method_name)(geometry, src_crs)
def _project_point(self, point, src_crs):
return sgeom.Point(*self.transform_point(point.x, point.y, src_crs))
def _project_line_string(self, geometry, src_crs):
return cartopy.trace.project_linear(geometry, src_crs, self)
def _project_linear_ring(self, linear_ring, src_crs):
"""
Project the given LinearRing from the src_crs into this CRS and
returns a list of LinearRings and a single MultiLineString.
"""
debug = False
# 1) Resolve the initial lines into projected segments
# 1abc
# def23ghi
# jkl41
multi_line_string = cartopy.trace.project_linear(linear_ring,
src_crs, self)
# Threshold for whether a point is close enough to be the same
# point as another.
threshold = max(np.abs(self.x_limits + self.y_limits)) * 1e-5
# 2) Simplify the segments where appropriate.
if len(multi_line_string.geoms) > 1:
# Stitch together segments which are close to continuous.
# This is important when:
# 1) The first source point projects into the map and the
# ring has been cut by the boundary.
# Continuing the example from above this gives:
# def23ghi
# jkl41abc
# 2) The cut ends of segments are too close to reliably
# place into an order along the boundary.
line_strings = list(multi_line_string.geoms)
any_modified = False
i = 0
if debug:
first_coord = np.array([ls.coords[0] for ls in line_strings])
last_coord = np.array([ls.coords[-1] for ls in line_strings])
print('Distance matrix:')
np.set_printoptions(precision=2)
x = first_coord[:, np.newaxis, :]
y = last_coord[np.newaxis, :, :]
print(np.abs(x - y).max(axis=-1))
while i < len(line_strings):
modified = False
j = 0
while j < len(line_strings):
if i != j and np.allclose(line_strings[i].coords[0],
line_strings[j].coords[-1],
atol=threshold):
if debug:
print(f'Joining together {i} and {j}.')
last_coords = list(line_strings[j].coords)
first_coords = list(line_strings[i].coords)[1:]
combo = sgeom.LineString(last_coords + first_coords)
if j < i:
i, j = j, i
del line_strings[j], line_strings[i]
line_strings.append(combo)
modified = True
any_modified = True
break
else:
j += 1
if not modified:
i += 1
if any_modified:
multi_line_string = sgeom.MultiLineString(line_strings)
# 3) Check for rings that have been created by the projection stage.
rings = []
line_strings = []
for line in multi_line_string.geoms:
if len(line.coords) > 3 and np.allclose(line.coords[0],
line.coords[-1],
atol=threshold):
result_geometry = sgeom.LinearRing(line.coords[:-1])
rings.append(result_geometry)
else:
line_strings.append(line)
# If we found any rings, then we should re-create the multi-line str.
if rings:
multi_line_string = sgeom.MultiLineString(line_strings)
return rings, multi_line_string
def _project_multipoint(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
geoms.append(self._project_point(geom, src_crs))
if geoms:
return sgeom.MultiPoint(geoms)
else:
return sgeom.MultiPoint()
def _project_multiline(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
r = self._project_line_string(geom, src_crs)
if r:
geoms.extend(r.geoms)
if geoms:
return sgeom.MultiLineString(geoms)
else:
return []
def _project_multipolygon(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
r = self._project_polygon(geom, src_crs)
if r:
geoms.extend(r.geoms)
if geoms:
result = sgeom.MultiPolygon(geoms)
else:
result = sgeom.MultiPolygon()
return result
def _project_polygon(self, polygon, src_crs):
"""
Return the projected polygon(s) derived from the given polygon.
"""
# Determine orientation of polygon.
# TODO: Consider checking the internal rings have the opposite
# orientation to the external rings?
if src_crs.is_geodetic():
is_ccw = True
else:
is_ccw = polygon.exterior.is_ccw
# Project the polygon exterior/interior rings.
# Each source ring will result in either a ring, or one or more
# lines.
rings = []
multi_lines = []
for src_ring in [polygon.exterior] + list(polygon.interiors):
p_rings, p_mline = self._project_linear_ring(src_ring, src_crs)
if p_rings:
rings.extend(p_rings)
if len(p_mline.geoms) > 0:
multi_lines.append(p_mline)
# Convert any lines to rings by attaching them to the boundary.
if multi_lines:
rings.extend(self._attach_lines_to_boundary(multi_lines, is_ccw))
# Resolve all the inside vs. outside rings, and convert to the
# final MultiPolygon.
return self._rings_to_multi_polygon(rings, is_ccw)
def _attach_lines_to_boundary(self, multi_line_strings, is_ccw):
"""
Return a list of LinearRings by attaching the ends of the given lines
to the boundary, paying attention to the traversal directions of the
lines and boundary.
"""
debug = False
debug_plot_edges = False
# Accumulate all the boundary and segment end points, along with
# their distance along the boundary.
edge_things = []
# Get the boundary as a LineString of the correct orientation
# so we can compute distances along it.
if is_ccw:
boundary = self.ccw_boundary
else:
boundary = self.cw_boundary