@@ -1061,32 +1061,39 @@ def __init__(
10611061 false_northing = None ,
10621062 true_scale_lat = None ,
10631063 ellipsoid = None ,
1064+ scale_factor_at_projection_origin = None ,
10641065 ):
10651066 """
10661067 Constructs a Stereographic coord system.
10671068
1068- Args:
1069+ Parameters
1070+ ----------
10691071
1070- * central_lat:
1072+ central_lat : float
10711073 The latitude of the pole.
10721074
1073- * central_lon:
1075+ central_lon : float
10741076 The central longitude, which aligns with the y axis.
10751077
1076- Kwargs:
1077-
1078- * false_easting:
1079- X offset from planar origin in metres. Defaults to 0.0 .
1078+ false_easting : float, optional
1079+ X offset from planar origin in metres.
10801080
1081- * false_northing:
1082- Y offset from planar origin in metres. Defaults to 0.0 .
1081+ false_northing : float, optional
1082+ Y offset from planar origin in metres.
10831083
1084- * true_scale_lat:
1084+ true_scale_lat : float, optional
10851085 Latitude of true scale.
10861086
1087- * ellipsoid (:class:`GeogCS`):
1087+ scale_factor_at_projection_origin : float, optional
1088+ Scale factor at the origin of the projection
1089+
1090+ ellipsoid : :class:`GeogCS`, optional
10881091 If given, defines the ellipsoid.
10891092
1093+ Notes
1094+ -----
1095+ It is only valid to provide one of true_scale_lat and scale_factor_at_projection_origin
1096+
10901097 """
10911098
10921099 #: True latitude of planar origin in degrees.
@@ -1105,27 +1112,42 @@ def __init__(
11051112 self .true_scale_lat = _arg_default (
11061113 true_scale_lat , None , cast_as = _float_or_None
11071114 )
1108- # N.B. the way we use this parameter, we need it to default to None,
1115+ #: Scale factor at projection origin.
1116+ self .scale_factor_at_projection_origin = _arg_default (
1117+ scale_factor_at_projection_origin , None , cast_as = _float_or_None
1118+ )
1119+ # N.B. the way we use these parameters, we need them to default to None,
11091120 # and *not* to 0.0 .
11101121
1122+ if (
1123+ self .true_scale_lat is not None
1124+ and self .scale_factor_at_projection_origin is not None
1125+ ):
1126+ raise ValueError (
1127+ "It does not make sense to provide both "
1128+ '"scale_factor_at_projection_origin" and "true_scale_latitude". '
1129+ )
1130+
11111131 #: Ellipsoid definition (:class:`GeogCS` or None).
11121132 self .ellipsoid = ellipsoid
11131133
1114- def __repr__ (self ):
1115- return (
1116- "Stereographic(central_lat={!r}, central_lon={!r}, "
1117- "false_easting={!r}, false_northing={!r}, "
1118- "true_scale_lat={!r}, "
1119- "ellipsoid={!r})" .format (
1120- self .central_lat ,
1121- self .central_lon ,
1122- self .false_easting ,
1123- self .false_northing ,
1124- self .true_scale_lat ,
1125- self .ellipsoid ,
1134+ def _repr_attributes (self ):
1135+ if self .scale_factor_at_projection_origin is None :
1136+ scale_info = "true_scale_lat={!r}, " .format (self .true_scale_lat )
1137+ else :
1138+ scale_info = "scale_factor_at_projection_origin={!r}, " .format (
1139+ self .scale_factor_at_projection_origin
11261140 )
1141+ return (
1142+ f"(central_lat={ self .central_lat } , central_lon={ self .central_lon } , "
1143+ f"false_easting={ self .false_easting } , false_northing={ self .false_northing } , "
1144+ f"{ scale_info } "
1145+ f"ellipsoid={ self .ellipsoid } )"
11271146 )
11281147
1148+ def __repr__ (self ):
1149+ return "Stereographic" + self ._repr_attributes ()
1150+
11291151 def as_cartopy_crs (self ):
11301152 globe = self ._ellipsoid_to_globe (self .ellipsoid , ccrs .Globe ())
11311153
@@ -1135,13 +1157,81 @@ def as_cartopy_crs(self):
11351157 self .false_easting ,
11361158 self .false_northing ,
11371159 self .true_scale_lat ,
1160+ self .scale_factor_at_projection_origin ,
11381161 globe = globe ,
11391162 )
11401163
11411164 def as_cartopy_projection (self ):
11421165 return self .as_cartopy_crs ()
11431166
11441167
1168+ class PolarStereographic (Stereographic ):
1169+ """
1170+ A subclass of the stereographic map projection centred on a pole.
1171+
1172+ """
1173+
1174+ grid_mapping_name = "polar_stereographic"
1175+
1176+ def __init__ (
1177+ self ,
1178+ central_lat ,
1179+ central_lon ,
1180+ false_easting = None ,
1181+ false_northing = None ,
1182+ true_scale_lat = None ,
1183+ scale_factor_at_projection_origin = None ,
1184+ ellipsoid = None ,
1185+ ):
1186+ """
1187+ Construct a Polar Stereographic coord system.
1188+
1189+ Parameters
1190+ ----------
1191+
1192+ central_lat : {90, -90}
1193+ The latitude of the pole.
1194+
1195+ central_lon : float
1196+ The central longitude, which aligns with the y axis.
1197+
1198+ false_easting : float, optional
1199+ X offset from planar origin in metres.
1200+
1201+ false_northing : float, optional
1202+ Y offset from planar origin in metres.
1203+
1204+ true_scale_lat : float, optional
1205+ Latitude of true scale.
1206+
1207+ scale_factor_at_projection_origin : float, optional
1208+ Scale factor at the origin of the projection
1209+
1210+ ellipsoid : :class:`GeogCS`, optional
1211+ If given, defines the ellipsoid.
1212+
1213+ Notes
1214+ -----
1215+ It is only valid to provide at most one of `true_scale_lat` and
1216+ `scale_factor_at_projection_origin`.
1217+
1218+
1219+ """
1220+
1221+ super ().__init__ (
1222+ central_lat = central_lat ,
1223+ central_lon = central_lon ,
1224+ false_easting = false_easting ,
1225+ false_northing = false_northing ,
1226+ true_scale_lat = true_scale_lat ,
1227+ scale_factor_at_projection_origin = scale_factor_at_projection_origin ,
1228+ ellipsoid = ellipsoid ,
1229+ )
1230+
1231+ def __repr__ (self ):
1232+ return "PolarStereographic" + self ._repr_attributes ()
1233+
1234+
11451235class LambertConformal (CoordSystem ):
11461236 """
11471237 A coordinate system in the Lambert Conformal conic projection.
0 commit comments