acgc.solar
Module for calculating solar position and TOA insolation
The functions here are are vectorized and generally broadcast over xarray dimensions, making this program faster than PySolar and pvlib. Calculations here use orbital parameters from the NOAA Solar Position Calculator, following Jean Meeus Astronomical Algorithms, unless the "fast" keyword is used. These accurate calculations are suitable for years -2000 to +3000. The "fast" calculations have lower accuracy orbital parameters and coarser approximation for the equation of time. All calculations here use geocentric solar position, neglecting the parallax effect of viewing the sun from different points on Earth (i.e. topocentric vs. geocentric in NREL SPA algorithm).
Accuracy: The NREL SPA algorithm in pvlib is used as an accurate reference. The maximum error in solar zenith angle computed here is 0.02° over 1900-2100. The maximum error in overall solar angular position is 0.022°. Large apparent differences in azimuth alone can occur when the sun is near zenith or nadir, where a small angular displacement results in a large azimuthal change.
The "fast" calculations have typical errors of ~0.2°.
NOAA Solar Calculator https://gml.noaa.gov/grad/solcalc/calcdetails.html
1#!/usr/local/bin/env python3 2'''Module for calculating solar position and TOA insolation 3 4The functions here are are vectorized and generally broadcast over xarray dimensions, 5making this program faster than PySolar and pvlib. Calculations here use orbital parameters 6from the NOAA Solar Position Calculator, following Jean Meeus Astronomical Algorithms, unless 7the "fast" keyword is used. These accurate calculations are suitable for years -2000 to +3000. 8The "fast" calculations have lower accuracy orbital parameters and coarser approximation 9for the equation of time. All calculations here use geocentric solar position, neglecting the 10parallax effect of viewing the sun from different points on Earth (i.e. topocentric vs. 11geocentric in NREL SPA algorithm). 12 13Accuracy: 14The NREL SPA algorithm in pvlib is used as an accurate reference. 15The maximum error in solar zenith angle computed here is 0.02° over 1900-2100. 16The maximum error in overall solar angular position is 0.022°. 17Large apparent differences in azimuth alone can occur when the sun is near zenith or nadir, 18where a small angular displacement results in a large azimuthal change. 19 20The "fast" calculations have typical errors of ~0.2°. 21 22NOAA Solar Calculator 23https://gml.noaa.gov/grad/solcalc/calcdetails.html 24''' 25 26from collections import namedtuple 27import warnings 28import numpy as np 29import pandas as pd 30import xarray as xr 31 32def insolation_toa( lat, lon, datetime, solar_pos=None, **kwargs ): 33 '''Insolation at top of the atmosphere, accounting for solar zenith angle 34 35 Parameters 36 ---------- 37 lat : float or ndarray 38 latitude in degrees 39 lon : float or ndarray 40 longitudes in degrees 41 datetime : datetime-like or str 42 date and time. Include time zone or UTC will be assumed 43 solar_pos : `SolarPositionResults`, optional 44 solar position parameters from a prior call to `solar_position`. 45 This can shorten runtime when calling multiple functions that use solar_position. 46 **kwargs passed to `solar_zenith_angle` 47 48 Returns 49 ------- 50 Insolation : float 51 radiation flux density accounting for solar zenith angle, W/m2 52 ''' 53 54 if solar_pos is None: 55 solar_pos = solar_position( datetime ) 56 57 S = solar_constant(datetime, solar_pos=solar_pos) 58 sza = solar_zenith_angle(lat, lon, datetime, **kwargs, solar_pos=solar_pos ) 59 60 return S * np.cos(sza) 61 62def solar_azimuth_angle( lat, lon, datetime, solar_pos=None ): 63 '''Solar azimuth angle (degrees) for a latitude, longitude, date and time 64 65 SAA is degrees clockwise from north. 66 67 Parameters 68 ---------- 69 lat : float or ndarray 70 latitude in degrees 71 lon : float or ndarray 72 longitudes in degrees 73 datetime : datetime-like or str 74 date and time. Include time zone or UTC will be assumed 75 solar_pos : `SolarPositionResults`, optional 76 solar position parameters from a prior call to `solar_position`. 77 This can shorten runtime when calling multiple functions that use solar_position. 78 79 Returns 80 ------- 81 saa : float or ndarray 82 solar azimuth angle in degrees (clockwise from north) 83 ''' 84 # Convert to pandas Timestamp, if needed 85 datetime = _to_datetime(datetime) 86 87 # Subsolar point, latitude longitude, degrees 88 solar_lat = solar_latitude( datetime, solar_pos=solar_pos ) 89 solar_lon = solar_longitude( datetime, solar_pos=solar_pos ) 90 91 # Vector pointing toward sun 92 x = np.cos( solar_lat * pi180 ) * np.sin( (solar_lon - lon) * pi180 ) 93 y = np.cos( lat*pi180 ) * np.sin( solar_lat*pi180 ) \ 94 - np.sin( lat*pi180 ) * np.cos( solar_lat*pi180 ) \ 95 * np.cos( (solar_lon - lon) * pi180 ) 96 97 # Azimuth angle from north, degrees 98 saa = np.arctan2( x, y ) / pi180 99 100 # Change range [-180,180] to [0,360] 101 return np.mod( saa+360, 360 ) 102 103def solar_elevation_angle( lat, lon, datetime, alt=0, 104 refraction=False, temperature=10., pressure=101325., 105 solar_pos=None ): 106 '''Solar elevation angle (degrees) above the horizon 107 108 The altitude parameter should be the vertical distance 109 above the surrounding terrain that defines the horizon, 110 not necessarily the altitude above sea level or the altitude above ground level. 111 For example, on a mountain peak that is 4000 m above sea level and 112 1500 m above the surrounding plateau, the relevant altitude is 1500 m. 113 For an observer on the plateau, the relevant altitude is 0 m. 114 115 See documentation for `solar_zenith_angle` and `horizon_zenith_angle`. 116 117 Parameters 118 ---------- 119 lat : float or ndarray 120 latitude in degrees 121 lon : float or ndarray 122 longitudes in degrees 123 datetime : datetime-like or str 124 date and time. Include time zone or UTC will be assumed 125 alt : float or ndarray (default=0) 126 altitude above surrounding terrain that defines the horizon, meters 127 refraction : bool (default=False) 128 specifies whether to account for atmospheric refraction 129 temperature : float or ndarray (default=10) 130 surface atmospheric temperature (Celsius), only used for refraction calculation 131 pressure : float or ndarray (default=101325) 132 surface atmospheric pressure (Pa), only used for refraction calculation 133 solar_pos : `SolarPositionResults`, optional 134 solar position parameters from a prior call to `solar_position`. 135 This can shorten runtime when calling multiple functions that use solar_position. 136 137 Returns 138 ------- 139 sea : float or ndarray 140 solar elevation angle in degrees at the designated locations and times 141 If refraction=False, this is the true solar elevation angle 142 If refraction=True, this is the apparent solar elevation angle 143 144 ''' 145 146 if refraction and np.any(alt): 147 warnings.warn( 'Atmospheric refraction is calculated for surface conditions, ' 148 + 'but an altitude above the surface was specified', 149 category=UserWarning, 150 stacklevel=2 ) 151 152 sea = horizon_zenith_angle( alt ) \ 153 - solar_zenith_angle( lat, lon, datetime, refraction, temperature, pressure, 154 solar_pos=solar_pos ) 155 156 return sea 157 158def solar_zenith_angle( lat, lon, datetime, 159 refraction=False, temperature=10., pressure=101325., 160 solar_pos=None ): 161 '''Solar zenith angle (degrees) for a given latitude, longitude, date and time. 162 163 Accounts for equation of time and (optionally) for atmospheric refraction. 164 Altitude of the observer is not accounted for, which can be important when the sun 165 is near the horizon. 166 167 Results are accurate to tenths of a degree, except where altitude is important 168 (< 20 degrees solar elevation) 169 170 Parameters 171 ---------- 172 lat : float or ndarray 173 latitude in degrees 174 lon : float or ndarray 175 longitudes in degrees 176 datetime : datetime-like or str 177 date and time. Include time zone or UTC will be assumed 178 refraction : bool (default=False) 179 specifies whether to account for atmospheric refraction 180 temperature : float or ndarray (default=10) 181 surface atmospheric temperature (Celsius), only used for refraction calculation 182 pressure : float or ndarray (default=101325) 183 surface atmospheric pressure (Pa), only used for refraction calculation 184 solar_pos : `SolarPositionResults`, optional 185 solar position parameters from a prior call to `solar_position`. 186 This can shorten runtime when calling multiple functions that use solar_position. 187 188 Returns 189 ------- 190 sza : float or ndarray 191 solar zenith angle in degrees at the designated locations and times 192 If refraction=False, this is the true solar zenith angle 193 If refraction=True, this is the apparent solar zenith angle 194 ''' 195 # Convert to pandas Timestamp, if needed 196 datetime = _to_datetime(datetime) 197 198 # Solar declination, degrees 199 if solar_pos is None: 200 dec = solar_declination( datetime ) 201 else: 202 dec = solar_pos.declination 203 204 # Hour angle, degrees 205 Ha = solar_hour_angle( lon, datetime, solar_pos ) 206 207 # True solar zenith angle, radians 208 sza = np.arccos( np.sin(lat*pi180) * np.sin(dec*pi180) + \ 209 np.cos(lat*pi180) * np.cos(dec*pi180) * np.cos(Ha*pi180) ) 210 211 # Convert radians -> degrees 212 sza /= pi180 213 214 if refraction: 215 # Subtract refraction angle (degrees) from zenith angle. 216 # SZA is always smaller due to refraction. 217 sza -= refraction_angle( 90-sza, pressure, temperature ) 218 219 return sza 220 221def sunrise_time( *args, **kwargs ): 222 '''Compute sunrise time 223 224 See `sun_times` for Parameters.''' 225 result = sun_times( *args, **kwargs ) 226 return result[0] 227 228def sunset_time( *args, **kwargs ): 229 '''Compute sunset time 230 231 See `sun_times` for Parameters.''' 232 result = sun_times( *args, **kwargs ) 233 return result[1] 234 235def day_length( *args, **kwargs ): 236 '''Compute length of daylight 237 238 See `sun_times` for Parameters.''' 239 result = sun_times( *args, **kwargs ) 240 return result[2] 241 242def solar_noon( *args, **kwargs ): 243 '''Compute time of solar noon (meridian transit) 244 245 See `sun_times` for Parameters.''' 246 result = sun_times( *args, **kwargs ) 247 return result[3] 248 249def solar_constant( datetime, solar_pos=None ): 250 '''Compute solar constant for specific date or dates 251 252 Parameters 253 ---------- 254 datetime : datetime-like 255 solar_pos : `SolarPositionResults`, optional 256 solar position parameters from a prior call to `solar_position`. 257 This can shorten runtime when calling multiple functions that use solar_position. 258 259 Returns 260 ------- 261 S : float 262 Solar direct beam radiation flux density, W/m2 263 ''' 264 if solar_pos is None: 265 solar_pos = solar_position( datetime ) 266 S = 1361/solar_pos.distance**2 267 268 return S 269 270def solar_declination( datetime, fast=False, solar_pos=None ): 271 '''Calculate solar declination (degrees) for specified date 272 273 Parameters 274 ---------- 275 datetime : datetime-like or str 276 date and time. Include time zone or UTC will be assumed 277 fast : bool (default=False) 278 Specifies using a faster but less accurate calculation 279 solar_pos : `SolarPositionResults`, optional 280 solar position parameters from a prior call to `solar_position`. 281 This can shorten runtime when calling multiple functions that use solar_position. 282 283 Returns 284 ------- 285 dec : float 286 solar declination in degrees at the specified date 287 ''' 288 # Convert to pandas Timestamp, if needed 289 datetime = _to_datetime(datetime) 290 291 # Select the accurate or fast calculation 292 accurate = not fast 293 294 if solar_pos is not None: 295 dec = solar_pos.declination 296 297 else: 298 if accurate: 299 300 # Solar declination, degrees 301 dec, junk, junk, junk, junk = solar_position( datetime ) 302 303 else: 304 # The fast method implements 305 # Eq. 9.68-9.72 from M.Z. Jacobson, Fundamentals of Atmospheric Modeling 306 307 # Number of days since beginning of 2000 308 NJD = datetime - np.datetime64('2000-01-01') 309 try: 310 NJD = NJD.dt.days 311 except AttributeError: 312 NJD = NJD.days 313 314 # Obliquity, degrees 315 ob = 23.439 - 4e-7 * NJD 316 317 # Parameters for ecliptic, degrees 318 gm = 357.528 + 0.9856003 * NJD 319 lm = 280.460 + 0.9856474 * NJD 320 321 # Ecliptic longitude of sun, degrees 322 ec = lm + 1.915 * np.sin( gm * pi180 ) + 0.020 * np.sin( 2 * gm * pi180 ) 323 324 #Solar declination, degrees 325 dec = np.arcsin( np.sin( ob * pi180 ) * np.sin( ec * pi180 ) ) / pi180 326 327 return dec 328 329def solar_latitude( *args, **kwargs ): 330 '''Latitude of the subsolar point, degrees 331 332 Parameters 333 ---------- 334 datetime : datetime-like or str 335 date and time. Include time zone or UTC will be assumed 336 fast : bool 337 see `solar_declination` 338 solar_pos : `SolarPositionResults`, optional 339 solar position parameters from a prior call to `solar_position`. 340 This can shorten runtime when calling multiple functions that use solar_position. 341 342 Returns 343 ------- 344 latitude : float 345 degrees of latitude 346 ''' 347 return solar_declination( *args, **kwargs ) 348 349def solar_longitude( datetime, solar_pos=None ): 350 '''Longitude of the subsolar point, degrees 351 352 Parameters 353 ---------- 354 datetimeUTC : datetime-like or str 355 date and time. Include time zone or UTC will be assumed 356 solar_pos : `SolarPositionResults`, optional 357 solar position parameters from a prior call to `solar_position`. 358 This can shorten runtime when calling multiple functions that use solar_position. 359 360 Returns 361 ------- 362 longitude : float 363 degrees of longitude 364 ''' 365 # Convert to pandas Timestamp, if needed 366 datetimeUTC, tz_in = _to_datetime_utc(datetime) 367 368 # Longitude of subsolar point, degrees 369 # Equation of time will be added below 370 try: 371 # Treat as xarray.DataArray or pandas.Series 372 solar_lon = - 15 * ( datetimeUTC.dt.hour + 373 datetimeUTC.dt.minute / 60 + 374 datetimeUTC.dt.second / 3600 - 12 ) 375 except AttributeError: 376 solar_lon = - 15 * ( datetimeUTC.hour + 377 datetimeUTC.minute / 60 + 378 datetimeUTC.second / 3600 - 12 ) 379 380 # Add equation of time to the solar longitude, degrees 381 if solar_pos is None: 382 eot = equation_of_time( datetimeUTC, degrees=True ) 383 else: 384 eot = solar_pos.equation_of_time 385 solar_lon -= eot 386 387 return solar_lon 388 389def solar_hour_angle( lon, datetime, solar_pos=None ): 390 '''Solar hour angle (degrees) for specified longitude, date and time 391 392 Hour angle is the angular displacement of the sun from the local meridian. 393 It is zero at local noon, negative in the morning, and positive is afternoon. 394 395 Parameters 396 ---------- 397 lon : float 398 longitude in degrees east 399 datetimeUTC : datetime-like or str 400 date and time. Include time zone or UTC will be assumed 401 solar_pos : `SolarPositionResults`, optional 402 solar position parameters from a prior call to `solar_position`. 403 This can shorten runtime when calling multiple functions that use solar_position. 404 405 Returns 406 ------- 407 ha : float 408 hour angle in degrees at the specified location and time 409 ''' 410 411 # Subsolar longitude, degrees 412 solar_lon = solar_longitude(datetime, solar_pos ) 413 414 # Hour angle, degrees 415 Ha = lon - solar_lon 416 417 return Ha 418 419def equation_of_time( datetime, degrees=False, fast=False ): 420 '''Equation of time for specified date 421 422 Accounts for the solar day being slightly different from 24 hours 423 424 Parameters 425 ---------- 426 datetime : datetime-like or str 427 date and time. Include time zone or UTC will be assumed 428 degrees : bool (default=False) 429 If True, then return value in compass degrees 430 If False, then return value in minutes of an hour 431 fast : bool (default=False) 432 specifies whether to use a faster, but less accurate calculation 433 434 Returns 435 ------- 436 eot : float 437 equation of time on the specified date, degrees or minutes 438 ''' 439 # Convert to pandas Timestamp, if needed 440 datetime = _to_datetime(datetime) 441 442 # Determine whether to use the fast or accurate calculation 443 accurate = not fast 444 445 if accurate: 446 447 # Equation of time, minutes 448 junk, junk, eot, junk, junk = solar_position( datetime ) 449 450 else: 451 # Implements the "alternative equation" from Wikipedia, derived from 452 # https://web.archive.org/web/20120323231813/http://www.green-life-innovators.org/tiki-index.php?page=The%2BLatitude%2Band%2BLongitude%2Bof%2Bthe%2BSun%2Bby%2BDavid%2BWilliams 453 # When compared to the NREL SPA algorithm, differences reach are up to about 0.5 minute. 454 # Note: Leap years are not accounted for. 455 456 # Equation of time, accounts for the solar day differing slightly from 24 hr 457 try: 458 doy = datetime.dt.dayofyear 459 except AttributeError: 460 doy = datetime.dayofyear 461 W = 360 / 365.24 462 A = W * (doy+10) 463 B = A + 1.914 * np.sin( W * (doy-2) * pi180 ) 464 C = ( A - np.arctan2( np.tan(B*pi180), np.cos(23.44*pi180) ) / pi180 ) / 180 465 466 # Equation of time in minutes of an hour (1440 minutes per day) 467 eot = 720 * ( C - np.round(C) ) 468 469 # Equation of time, minutes -> degrees (360 degrees per day) 470 if degrees: 471 eot = eot / 60 * 360 / 24 472 473 return eot 474 475def solar_position( datetime ): 476 '''Compute position of sun (declination, right ascension, equation of time, distance)) on specified date 477 478 Calculations follow the NOAA solar calculator spreadsheet 479 Applicable to years 1900-2100. 480 481 Parameters 482 ---------- 483 date : datetime-like or str 484 date and time. Include time zone or UTC will be assumed 485 486 Returns 487 ------- 488 result : `SolarPositionResults` 489 declination, right ascension, equation of time, Earth-sun distance, datetime 490 ''' 491 # Ensure time is Timestamp in UTC 492 datetimeUTC, tz_in = _to_datetime_utc(datetime) 493 494 # Raise warning if any dates are outside date range 495 # recommended for orbital parameters used here 496 if np.logical_or( np.any( datetimeUTC < np.datetime64('1900-01-01') ), 497 np.any( datetimeUTC > np.datetime64('2100-01-01') ) ): 498 warnings.warn('Solar position accuracy declines for dates outside 1900-2100', \ 499 RuntimeWarning ) 500 501 # Number of days since 1 Jan 2000 502 NJD = datetimeUTC - np.datetime64('2000-01-01') 503 try: 504 NJD = NJD.dt.days \ 505 + NJD.dt.seconds / 86400. 506 except AttributeError: 507 NJD = NJD.days \ 508 + NJD.seconds / 86400. 509 510 # Julian day (since 12:00UTC 1 Jan 4713 BCE) 511 NJD += 2451544.50 512 513 # Julian century 514 JC = (NJD-2451545)/36525 515 516 # Earth orbital eccentricity, unitless 517 ec = 0.016708634 - JC*( 0.000042037 + 0.0000001267*JC ) 518 519 # Earth mean orbital obliquity, degree 520 mean_ob = 23 + ( 26 + ( (21.448 521 - JC * (46.815 522 + JC * (0.00059 - JC * 0.001813) ) ) )/60 )/60 523 524 # Earth true orbital obliquity, corrected for nutation, degree 525 ob = mean_ob + 0.00256 * np.cos( (125.04 - 1934.136*JC ) * pi180 ) 526 527 # Sun Mean ecliptic longitude, degree 528 mean_ec_lon = np.mod( 280.46646 + JC*( 36000.76983 + JC*0.0003032 ), 360 ) 529 530 # Sun Mean anomaly, degree 531 mean_anom = 357.52911 + JC*( 35999.05029 - 0.0001537*JC ) 532 533 # Sun Equation of center, degree 534 eq_center = np.sin(mean_anom*pi180) * (1.914602 - JC*( 0.004817 + 0.000014*JC )) \ 535 + np.sin(2*mean_anom*pi180) * (0.019993 - 0.000101*JC) \ 536 + np.sin(3*mean_anom*pi180) * 0.000289 537 538 # Sun True ecliptic longitude, degrees 539 true_ec_lon = mean_ec_lon + eq_center 540 541 # Sun True anomaly, degree 542 true_anom = mean_anom + eq_center 543 544 # Earth-Sun distance, AU 545 distance = (1.000001018 * (1-ec**2) ) / (1 + ec * np.cos( true_anom * pi180 )) 546 547 # Sun Apparent ecliptic longitude, corrected for nutation, degrees 548 ec_lon = true_ec_lon - 0.00569 - 0.00478 * np.sin( (125.04 - 1934.136*JC ) * pi180) 549 550 # Sun Right ascension, deg 551 right_ascension = np.arctan2( np.cos(ob*pi180) * np.sin(ec_lon*pi180), 552 np.cos(ec_lon*pi180) ) / pi180 553 554 # Sun Declination, deg 555 declination = np.arcsin( np.sin(ob*pi180) * np.sin(ec_lon*pi180) ) / pi180 556 557 # var y 558 vary = np.tan( ob/2 * pi180 )**2 559 560 # Equation of time, minutes 561 eot = vary * np.sin( 2 * mean_ec_lon * pi180) \ 562 - 2 * ec * np.sin( mean_anom * pi180 ) \ 563 + 4 * ec * vary * np.sin( mean_anom * pi180 ) * np.cos( 2 * mean_ec_lon * pi180) \ 564 - 0.5 * vary**2 * np.sin( 4 * mean_ec_lon * pi180) \ 565 - 1.25 * ec**2 * np.sin( 2 * mean_anom * pi180) 566 eot = eot * 4 / pi180 567 568 # Bundle results 569 result = SolarPositionResults(declination, 570 right_ascension, 571 eot, 572 distance, 573 datetimeUTC) 574 575 return result 576 577def sun_times( lat, lon, datetime, tz_out=None, sza_sunrise=90.833, 578 fast=False, solar_pos=None ): 579 '''Compute times of sunrise, sunset, solar noon, and day length 580 581 Common options for solar zenith angle at sunrise 582 - 90.833 for first edge of sun rising, typical (0.567°) refraction (default) 583 - 90.267 for first edge of sun rising, no refraction 584 - 90 degrees for center of sun rising, no refraction 585 586 Note: horizon_zenith_angle can be used to compute a more accurate horizon location 587 for sites at elevation. 588 589 Parameters 590 ---------- 591 lat : float or ndarray 592 latitude in degrees 593 lon : float or ndarray 594 longitudes in degrees 595 datetime : datetime-like or str 596 datetime, provide a time zone or UTC will be assumed 597 tz_out : str, pytz.timezone, datetime.tzinfo, optional 598 timezone to be used for output times. 599 If None is provided, then result will be in same time zone as input or UTC 600 sza_sunrise : float (default=90.833) 601 Solar zenith angle at which sunrise and sunset are calculated, degrees 602 fast : bool (default=False) 603 Select a faster but less accurate calculation 604 solar_pos : `SolarPositionResults`, optional 605 solar position parameters from a prior call to `solar_position`. 606 This can shorten runtime when calling multiple functions that use solar_position. 607 608 Returns 609 ------- 610 result : `SunTimesResults` 611 times of sunrise, sunset, solar noon, and day length 612 ''' 613 # Convert to pandas Timestamp in UTC, if needed 614 datetimeUTC, tz_in = _to_datetime_utc(datetime) 615 616 # If no output timezone is specified, use the input time zone 617 if tz_out is None: 618 tz_out = tz_in 619 620 # Select fast or accurate calculation 621 accurate = not fast 622 623 # Solar declination (degrees) and equation of time (minutes) 624 if solar_pos is not None: 625 dec = solar_pos.declination 626 eot = solar_pos.eot 627 else: 628 if accurate: 629 dec, junk, eot, junk, junk = solar_position( datetimeUTC ) 630 else: 631 dec = solar_declination( datetimeUTC ) 632 eot = equation_of_time( datetimeUTC ) 633 634 # Sunrise hour angle, degree 635 # Degrees east of the local meridian where sun rises 636 ha_sunrise = np.arccos( np.cos(sza_sunrise*pi180) / 637 (np.cos(lat*pi180) * np.cos(dec*pi180)) 638 - np.tan(lat*pi180)*np.tan(dec*pi180) ) / pi180 639 640 # Solar noon, local standard time, day fraction 641 solar_noon = (720 - 4*lon - eot ) / 1440 642 643 # Sunrise and sunset, local standard time, day fraction 644 t_sunrise = solar_noon - 4 * ha_sunrise / 1440 645 t_sunset = solar_noon + 4 * ha_sunrise / 1440 646 647 # Midnight UTC 648 # datetimeUTC is in UTC but time-zone-naive 649 try: 650 # Series time objects 651 dateUTC = datetimeUTC.dt.floor('D') 652 except AttributeError: 653 # Scalar time objects 654 dateUTC = datetimeUTC.floor('D') 655 656 # Convert day fraction -> date time 657 solar_noon = dateUTC + solar_noon * pd.Timedelta( 1, 'day' ) 658 t_sunrise = dateUTC + t_sunrise * pd.Timedelta( 1, 'day' ) 659 t_sunset = dateUTC + t_sunset * pd.Timedelta( 1, 'day' ) 660 661 # Convert to output timezone, if any is provided 662 if tz_out is not None: 663 if isinstance(solar_noon,(xr.DataArray,np.ndarray)) or \ 664 isinstance(t_sunrise,(xr.DataArray,np.ndarray)): 665 # These types don't localize tz, but we can add offset to the tz-naive time 666 if hasattr(datetimeUTC,'tz_localize'): 667 # For scalar datetime, there is a single time offset, which we can add 668 utcoffset = np.timedelta64( datetimeUTC.tz_localize(tz_out).utcoffset() ) 669 solar_noon += utcoffset 670 t_sunrise += utcoffset 671 t_sunset += utcoffset 672 else: 673 # For Series datetime, there are potentially multiple offsets. We can only add one 674 unique_datetimeUTC = pd.DatetimeIndex(np.unique(datetimeUTC)) 675 unique_utcoffsets = np.unique( unique_datetimeUTC.tz_localize('UTC') \ 676 - unique_datetimeUTC.tz_localize(tz_out) ) 677 if len(unique_utcoffsets)==1: 678 utcoffset = unique_utcoffsets[0] 679 solar_noon += utcoffset 680 t_sunrise += utcoffset 681 t_sunset += utcoffset 682 else: 683 # We might be able to handle multiple offsets if we can 684 # determine which dimension of `solar_noon` and `t_sunrise` 685 # is the time dimension. 686 raise ValueError('Multiple timezone offsets not supported. ' 687 +'Request output in UTC or reduce number of input times.') 688 else: 689 try: 690 # Series time objects 691 solar_noon = solar_noon.dt.tz_localize('UTC')\ 692 .dt.tz_convert(tz_out) 693 t_sunrise = t_sunrise.dt.tz_localize('UTC')\ 694 .dt.tz_convert(tz_out) 695 t_sunset = t_sunset.dt.tz_localize('UTC')\ 696 .dt.tz_convert(tz_out) 697 except AttributeError: 698 # Scale time objects 699 solar_noon = solar_noon.tz_localize('UTC')\ 700 .tz_convert(tz_out) 701 t_sunrise = t_sunrise.tz_localize('UTC')\ 702 .tz_convert(tz_out) 703 t_sunset = t_sunset.tz_localize('UTC')\ 704 .tz_convert(tz_out) 705 706 # Sunlight duration, minutes 707 day_length = 8 * ha_sunrise * pd.Timedelta(1, 'minute') 708 709 result = SunTimesResults(t_sunrise, t_sunset, day_length, solar_noon, datetimeUTC) 710 711 return result 712 713def horizon_zenith_angle( alt ): 714 '''Angle from the zenith to the horizon 715 716 The horizon is the locii of points where a line from the 717 observation location to the ellipsoid is tangent to the ellipsoid surface. 718 719 The altitude parameter should be the vertical distance 720 above the surrounding terrain that defines the horizon, 721 not necessarily the altitude above sea level or the altitude above ground level. 722 For example, on a mountain peak that is 4000 m above sea level and 723 1500 m above the surrounding plateau, the relevant altitude is 1500 m. 724 For an observer on the plateau, the relevant altitude is 0 m. 725 726 The implementation below assumes a spherical Earth. 727 Results using the WGS84 ellipsoid (see commented code below) 728 differ from the spherical case by << 1°. Terrain, 729 which is neglected here, has a larger effect on the horizon 730 location, so the simpler spherical calculation is appropriate. 731 732 Parameters 733 ---------- 734 lat : float or ndarray 735 latitude in degrees 736 alt : float or ndarray 737 altitude above surrounding terrain that defines the horizon, meters 738 739 Returns 740 ------- 741 hza : float or ndarray 742 horizon zenith angle in degrees 743 ''' 744 745 # WGS84 ellipsoid parameters 746 # semi-major radius, m 747 r_earth = 6378137.0 748 # ellipsoidal flattening, unitless 749 f = 1/298.257223563 750 751 # Horizon zenith angle, degrees (spherical earth) 752 hza = 180 - np.arcsin( r_earth / ( r_earth + alt ) ) / pi180 753 754 ## Ellipsoidal Earth 755 # # Eccentricity of ellipsoid 756 # ecc = f * (2-f) 757 # # Local (i.e. prime vertical) radius of curvature at latitude 758 # N = r_earth / np.sqrt( 1 - ecc**2 * np.sin(lat*pi180)**2 ) 759 # # Horizon zenith angle, degrees 760 # hza = 180 - np.arcsin( N / (N+alt) ) / pi180 761 762 return hza 763 764def refraction_angle( true_elevation_angle, pressure=101325., temperature_celsius=10. ): 765 '''Atmospheric refraction angle for light passing through Earth's atmosphere 766 767 The apparent locations in the sky of objects outsides Earth's atmosphere 768 differs from their true locations due to atmospheric refraction. 769 (e.g. sun and moon when they rise and set) 770 The apparent elevation of an object above the horizon is 771 apparent elevation angle = (true elevation angle) + (refraction angle) 772 773 The equations here are from Saemundsson/Bennett, whose calculations use 774 a typical vertical profile of atmospheric density (i.e. temperature and pressure). 775 The profiles can be rescaled to a particular surface temperature and pressure 776 to approximately account for varying atmospheric conditions. 777 Accurate refraction calculations should use fully specified vertical profile 778 of temperature and pressure, which cannot be done here. 779 780 Parameters 781 ---------- 782 true_elevation_angle : float 783 degrees above horizon of sun or other object 784 pressure : float (default=101325) 785 surface atmospheric pressure (Pa) 786 temperature_celsius : float (default=10) 787 surface atmospheric temperature (C) 788 789 Returns 790 ------- 791 angle : float 792 refraction angle in degrees. Value is zero when apparent elevation is below horizon 793 ''' 794 # Refraction angle, arcminutes 795 R = 1.02 / np.tan( ( true_elevation_angle + 10.3 / (true_elevation_angle + 5.11) ) * pi180 ) 796 # Account for temperature and pressure, arcminutes 797 R = R * pressure / 101325 * 283 / ( 273 + temperature_celsius ) 798 # Convert arcminutes -> degrees 799 R /= 60 800 801 # Result must be positive 802 R = np.maximum(R,0) 803 804 # Refraction defined only when the apparent elevation angle is positive 805 # Set refraction to zero when the apparent elevation is below horizon 806 refraction_angle = np.where( true_elevation_angle + R <= 0, 0, R) 807 808 return refraction_angle 809 810def _to_datetime(time_in): 811 '''Convert to Timestamp, Series of datetime64, or DataArray of datetime64 812 813 Parameters 814 ---------- 815 time_in : datetime-like or array 816 time to be converted 817 818 Returns 819 ------- 820 time_out : pandas.Timestamp, pandas.Series of datetime64, DataArray of datetime64 821 ''' 822 if hasattr(time_in,'dt'): 823 time_out = time_in 824 # elif isinstance(time_in, pd.DatetimeIndex ): 825 # # Unnecessary; DatetimeIndex will work fine in the else cases 826 # time_out = pd.Series(time_in) 827 # tz = time_out.dt.tz 828 else: 829 try: 830 # Convert list of times 831 time_out = pd.Series(pd.DatetimeIndex(time_in)) 832 except TypeError: 833 # Single datetime or str 834 time_out = pd.Timestamp(time_in) 835 836 return time_out 837 838def _to_datetime_utc( datetime_in ): 839 '''Convert to Timestamp, Series of datetime64, or DataArray of datetime64 and in UTC 840 841 Parameters 842 ---------- 843 datetime_in : datetime-like 844 date and time to be converted 845 846 Returns 847 ------- 848 datetimeUTC : Timestamp, Series, or DataArray 849 date and time in UTC but tz-naive 850 tz_in : datetime.timezone 851 timezone of datetime_in 852 ''' 853 854 # Ensure input is a timestamp 855 datetime_in = _to_datetime( datetime_in ) 856 857 # Convert to UTC, then strip timezone 858 try: 859 if hasattr(datetime_in,'dt'): 860 # Pandas Series objects 861 tz_in = datetime_in.dt.tz 862 datetimeUTC = datetime_in.dt.tz_convert('UTC').dt.tz_localize(None) 863 else: 864 # Scalar time objects 865 tz_in = datetime_in.tzinfo 866 datetimeUTC = datetime_in.tz_convert('UTC').tz_localize(None) 867 except (TypeError, AttributeError): 868 # tz-naive time objects: Timestamp, Series, (all) DataArrays 869 # No timezone info, so assume it is already UTC 870 warnings.warn('Time does not contain timezone. UTC will be assumed',RuntimeWarning) 871 datetimeUTC = datetime_in 872 tz_in = None 873 874 return datetimeUTC, tz_in 875 876SolarPositionResults = namedtuple('SolarPositionResults', 877 'declination right_ascension equation_of_time distance datetimeUTC') 878'''Namedtuple containing results of `solar_position` 879 880All values are geocentric, strictly accurate for the center of the Earth, not a point on 881Earth's surface. The parallax angle from Earth's center to surface is 4e-5 degrees. 882 883Attributes 884---------- 885declination : float 886 position of the sun relative to Earth's equatorial plane, degrees 887right_ascension : float 888 position of the sun along Earth's equatorial plane relative to the vernal equinox, degrees 889equation_of_time : float 890 equation of time (minutes) between mean solar time and true solar time. 891 Divide by 4 minutes per degree to obtain equation of time in degrees. 892distance : float 893 Earth-sun distance in AU (1 AU = 1.495978707e11 m) 894datetimeUTC : Timestamp, Series, or DataArray of numpy.datetime64 895 date and time input for calculations, UTC 896''' 897 898SunTimesResults = namedtuple('SunTimesResults', 899 'sunrise sunset day_length solar_noon datetimeUTC') 900'''Namedtuple containing results of `sun_times` 901 902Attributes 903---------- 904sunrise : pandas.DatetimeIndex 905 sunrise time, UTC if not specified otherwise 906sunset : pandas.DatetimeIndex 907 sunset time, UTC if not specified otherwise 908day_length : pandas.Timedelta 909 duration of daylight 910solar_noon : pandas.DatetimeIndex 911 time of meridian transit, UTC if not specified otherwise 912datetimeUTC : Timestamp, Series, or DataArray of numpy.datetime64 913 date and time input for calculations, UTC 914''' 915 916pi180 = np.pi / 180 917'''Constant $\\pi/180$''' 918 919# Aliases for functions 920sza = solar_zenith_angle 921'''Alias for `solar_zenith_angle`''' 922saa = solar_azimuth_angle 923'''Alias for `solar_azimuth_angle`''' 924sea = solar_elevation_angle 925'''Alias for `solar_elevation_angle`''' 926 927# Deprecated aliases 928equationOfTime = equation_of_time 929'''Alias for `equation_of_time`''' 930solarDeclination = solar_declination 931'''Alias for `solar_declination`'''
33def insolation_toa( lat, lon, datetime, solar_pos=None, **kwargs ): 34 '''Insolation at top of the atmosphere, accounting for solar zenith angle 35 36 Parameters 37 ---------- 38 lat : float or ndarray 39 latitude in degrees 40 lon : float or ndarray 41 longitudes in degrees 42 datetime : datetime-like or str 43 date and time. Include time zone or UTC will be assumed 44 solar_pos : `SolarPositionResults`, optional 45 solar position parameters from a prior call to `solar_position`. 46 This can shorten runtime when calling multiple functions that use solar_position. 47 **kwargs passed to `solar_zenith_angle` 48 49 Returns 50 ------- 51 Insolation : float 52 radiation flux density accounting for solar zenith angle, W/m2 53 ''' 54 55 if solar_pos is None: 56 solar_pos = solar_position( datetime ) 57 58 S = solar_constant(datetime, solar_pos=solar_pos) 59 sza = solar_zenith_angle(lat, lon, datetime, **kwargs, solar_pos=solar_pos ) 60 61 return S * np.cos(sza)
Insolation at top of the atmosphere, accounting for solar zenith angle
Parameters
- lat (float or ndarray): latitude in degrees
- lon (float or ndarray): longitudes in degrees
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position. - **kwargs passed to
solar_zenith_angle
Returns
- Insolation (float): radiation flux density accounting for solar zenith angle, W/m2
63def solar_azimuth_angle( lat, lon, datetime, solar_pos=None ): 64 '''Solar azimuth angle (degrees) for a latitude, longitude, date and time 65 66 SAA is degrees clockwise from north. 67 68 Parameters 69 ---------- 70 lat : float or ndarray 71 latitude in degrees 72 lon : float or ndarray 73 longitudes in degrees 74 datetime : datetime-like or str 75 date and time. Include time zone or UTC will be assumed 76 solar_pos : `SolarPositionResults`, optional 77 solar position parameters from a prior call to `solar_position`. 78 This can shorten runtime when calling multiple functions that use solar_position. 79 80 Returns 81 ------- 82 saa : float or ndarray 83 solar azimuth angle in degrees (clockwise from north) 84 ''' 85 # Convert to pandas Timestamp, if needed 86 datetime = _to_datetime(datetime) 87 88 # Subsolar point, latitude longitude, degrees 89 solar_lat = solar_latitude( datetime, solar_pos=solar_pos ) 90 solar_lon = solar_longitude( datetime, solar_pos=solar_pos ) 91 92 # Vector pointing toward sun 93 x = np.cos( solar_lat * pi180 ) * np.sin( (solar_lon - lon) * pi180 ) 94 y = np.cos( lat*pi180 ) * np.sin( solar_lat*pi180 ) \ 95 - np.sin( lat*pi180 ) * np.cos( solar_lat*pi180 ) \ 96 * np.cos( (solar_lon - lon) * pi180 ) 97 98 # Azimuth angle from north, degrees 99 saa = np.arctan2( x, y ) / pi180 100 101 # Change range [-180,180] to [0,360] 102 return np.mod( saa+360, 360 )
Solar azimuth angle (degrees) for a latitude, longitude, date and time
SAA is degrees clockwise from north.
Parameters
- lat (float or ndarray): latitude in degrees
- lon (float or ndarray): longitudes in degrees
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- saa (float or ndarray): solar azimuth angle in degrees (clockwise from north)
104def solar_elevation_angle( lat, lon, datetime, alt=0, 105 refraction=False, temperature=10., pressure=101325., 106 solar_pos=None ): 107 '''Solar elevation angle (degrees) above the horizon 108 109 The altitude parameter should be the vertical distance 110 above the surrounding terrain that defines the horizon, 111 not necessarily the altitude above sea level or the altitude above ground level. 112 For example, on a mountain peak that is 4000 m above sea level and 113 1500 m above the surrounding plateau, the relevant altitude is 1500 m. 114 For an observer on the plateau, the relevant altitude is 0 m. 115 116 See documentation for `solar_zenith_angle` and `horizon_zenith_angle`. 117 118 Parameters 119 ---------- 120 lat : float or ndarray 121 latitude in degrees 122 lon : float or ndarray 123 longitudes in degrees 124 datetime : datetime-like or str 125 date and time. Include time zone or UTC will be assumed 126 alt : float or ndarray (default=0) 127 altitude above surrounding terrain that defines the horizon, meters 128 refraction : bool (default=False) 129 specifies whether to account for atmospheric refraction 130 temperature : float or ndarray (default=10) 131 surface atmospheric temperature (Celsius), only used for refraction calculation 132 pressure : float or ndarray (default=101325) 133 surface atmospheric pressure (Pa), only used for refraction calculation 134 solar_pos : `SolarPositionResults`, optional 135 solar position parameters from a prior call to `solar_position`. 136 This can shorten runtime when calling multiple functions that use solar_position. 137 138 Returns 139 ------- 140 sea : float or ndarray 141 solar elevation angle in degrees at the designated locations and times 142 If refraction=False, this is the true solar elevation angle 143 If refraction=True, this is the apparent solar elevation angle 144 145 ''' 146 147 if refraction and np.any(alt): 148 warnings.warn( 'Atmospheric refraction is calculated for surface conditions, ' 149 + 'but an altitude above the surface was specified', 150 category=UserWarning, 151 stacklevel=2 ) 152 153 sea = horizon_zenith_angle( alt ) \ 154 - solar_zenith_angle( lat, lon, datetime, refraction, temperature, pressure, 155 solar_pos=solar_pos ) 156 157 return sea
Solar elevation angle (degrees) above the horizon
The altitude parameter should be the vertical distance above the surrounding terrain that defines the horizon, not necessarily the altitude above sea level or the altitude above ground level. For example, on a mountain peak that is 4000 m above sea level and 1500 m above the surrounding plateau, the relevant altitude is 1500 m. For an observer on the plateau, the relevant altitude is 0 m.
See documentation for solar_zenith_angle
and horizon_zenith_angle
.
Parameters
- lat (float or ndarray): latitude in degrees
- lon (float or ndarray): longitudes in degrees
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- alt (float or ndarray (default=0)): altitude above surrounding terrain that defines the horizon, meters
- refraction (bool (default=False)): specifies whether to account for atmospheric refraction
- temperature (float or ndarray (default=10)): surface atmospheric temperature (Celsius), only used for refraction calculation
- pressure (float or ndarray (default=101325)): surface atmospheric pressure (Pa), only used for refraction calculation
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- sea (float or ndarray): solar elevation angle in degrees at the designated locations and times If refraction=False, this is the true solar elevation angle If refraction=True, this is the apparent solar elevation angle
159def solar_zenith_angle( lat, lon, datetime, 160 refraction=False, temperature=10., pressure=101325., 161 solar_pos=None ): 162 '''Solar zenith angle (degrees) for a given latitude, longitude, date and time. 163 164 Accounts for equation of time and (optionally) for atmospheric refraction. 165 Altitude of the observer is not accounted for, which can be important when the sun 166 is near the horizon. 167 168 Results are accurate to tenths of a degree, except where altitude is important 169 (< 20 degrees solar elevation) 170 171 Parameters 172 ---------- 173 lat : float or ndarray 174 latitude in degrees 175 lon : float or ndarray 176 longitudes in degrees 177 datetime : datetime-like or str 178 date and time. Include time zone or UTC will be assumed 179 refraction : bool (default=False) 180 specifies whether to account for atmospheric refraction 181 temperature : float or ndarray (default=10) 182 surface atmospheric temperature (Celsius), only used for refraction calculation 183 pressure : float or ndarray (default=101325) 184 surface atmospheric pressure (Pa), only used for refraction calculation 185 solar_pos : `SolarPositionResults`, optional 186 solar position parameters from a prior call to `solar_position`. 187 This can shorten runtime when calling multiple functions that use solar_position. 188 189 Returns 190 ------- 191 sza : float or ndarray 192 solar zenith angle in degrees at the designated locations and times 193 If refraction=False, this is the true solar zenith angle 194 If refraction=True, this is the apparent solar zenith angle 195 ''' 196 # Convert to pandas Timestamp, if needed 197 datetime = _to_datetime(datetime) 198 199 # Solar declination, degrees 200 if solar_pos is None: 201 dec = solar_declination( datetime ) 202 else: 203 dec = solar_pos.declination 204 205 # Hour angle, degrees 206 Ha = solar_hour_angle( lon, datetime, solar_pos ) 207 208 # True solar zenith angle, radians 209 sza = np.arccos( np.sin(lat*pi180) * np.sin(dec*pi180) + \ 210 np.cos(lat*pi180) * np.cos(dec*pi180) * np.cos(Ha*pi180) ) 211 212 # Convert radians -> degrees 213 sza /= pi180 214 215 if refraction: 216 # Subtract refraction angle (degrees) from zenith angle. 217 # SZA is always smaller due to refraction. 218 sza -= refraction_angle( 90-sza, pressure, temperature ) 219 220 return sza
Solar zenith angle (degrees) for a given latitude, longitude, date and time.
Accounts for equation of time and (optionally) for atmospheric refraction. Altitude of the observer is not accounted for, which can be important when the sun is near the horizon.
Results are accurate to tenths of a degree, except where altitude is important (< 20 degrees solar elevation)
Parameters
- lat (float or ndarray): latitude in degrees
- lon (float or ndarray): longitudes in degrees
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- refraction (bool (default=False)): specifies whether to account for atmospheric refraction
- temperature (float or ndarray (default=10)): surface atmospheric temperature (Celsius), only used for refraction calculation
- pressure (float or ndarray (default=101325)): surface atmospheric pressure (Pa), only used for refraction calculation
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- sza (float or ndarray): solar zenith angle in degrees at the designated locations and times If refraction=False, this is the true solar zenith angle If refraction=True, this is the apparent solar zenith angle
222def sunrise_time( *args, **kwargs ): 223 '''Compute sunrise time 224 225 See `sun_times` for Parameters.''' 226 result = sun_times( *args, **kwargs ) 227 return result[0]
Compute sunrise time
See sun_times
for Parameters.
229def sunset_time( *args, **kwargs ): 230 '''Compute sunset time 231 232 See `sun_times` for Parameters.''' 233 result = sun_times( *args, **kwargs ) 234 return result[1]
Compute sunset time
See sun_times
for Parameters.
236def day_length( *args, **kwargs ): 237 '''Compute length of daylight 238 239 See `sun_times` for Parameters.''' 240 result = sun_times( *args, **kwargs ) 241 return result[2]
Compute length of daylight
See sun_times
for Parameters.
243def solar_noon( *args, **kwargs ): 244 '''Compute time of solar noon (meridian transit) 245 246 See `sun_times` for Parameters.''' 247 result = sun_times( *args, **kwargs ) 248 return result[3]
Compute time of solar noon (meridian transit)
See sun_times
for Parameters.
250def solar_constant( datetime, solar_pos=None ): 251 '''Compute solar constant for specific date or dates 252 253 Parameters 254 ---------- 255 datetime : datetime-like 256 solar_pos : `SolarPositionResults`, optional 257 solar position parameters from a prior call to `solar_position`. 258 This can shorten runtime when calling multiple functions that use solar_position. 259 260 Returns 261 ------- 262 S : float 263 Solar direct beam radiation flux density, W/m2 264 ''' 265 if solar_pos is None: 266 solar_pos = solar_position( datetime ) 267 S = 1361/solar_pos.distance**2 268 269 return S
Compute solar constant for specific date or dates
Parameters
datetime (datetime-like):
solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- S (float): Solar direct beam radiation flux density, W/m2
271def solar_declination( datetime, fast=False, solar_pos=None ): 272 '''Calculate solar declination (degrees) for specified date 273 274 Parameters 275 ---------- 276 datetime : datetime-like or str 277 date and time. Include time zone or UTC will be assumed 278 fast : bool (default=False) 279 Specifies using a faster but less accurate calculation 280 solar_pos : `SolarPositionResults`, optional 281 solar position parameters from a prior call to `solar_position`. 282 This can shorten runtime when calling multiple functions that use solar_position. 283 284 Returns 285 ------- 286 dec : float 287 solar declination in degrees at the specified date 288 ''' 289 # Convert to pandas Timestamp, if needed 290 datetime = _to_datetime(datetime) 291 292 # Select the accurate or fast calculation 293 accurate = not fast 294 295 if solar_pos is not None: 296 dec = solar_pos.declination 297 298 else: 299 if accurate: 300 301 # Solar declination, degrees 302 dec, junk, junk, junk, junk = solar_position( datetime ) 303 304 else: 305 # The fast method implements 306 # Eq. 9.68-9.72 from M.Z. Jacobson, Fundamentals of Atmospheric Modeling 307 308 # Number of days since beginning of 2000 309 NJD = datetime - np.datetime64('2000-01-01') 310 try: 311 NJD = NJD.dt.days 312 except AttributeError: 313 NJD = NJD.days 314 315 # Obliquity, degrees 316 ob = 23.439 - 4e-7 * NJD 317 318 # Parameters for ecliptic, degrees 319 gm = 357.528 + 0.9856003 * NJD 320 lm = 280.460 + 0.9856474 * NJD 321 322 # Ecliptic longitude of sun, degrees 323 ec = lm + 1.915 * np.sin( gm * pi180 ) + 0.020 * np.sin( 2 * gm * pi180 ) 324 325 #Solar declination, degrees 326 dec = np.arcsin( np.sin( ob * pi180 ) * np.sin( ec * pi180 ) ) / pi180 327 328 return dec
Calculate solar declination (degrees) for specified date
Parameters
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- fast (bool (default=False)): Specifies using a faster but less accurate calculation
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- dec (float): solar declination in degrees at the specified date
330def solar_latitude( *args, **kwargs ): 331 '''Latitude of the subsolar point, degrees 332 333 Parameters 334 ---------- 335 datetime : datetime-like or str 336 date and time. Include time zone or UTC will be assumed 337 fast : bool 338 see `solar_declination` 339 solar_pos : `SolarPositionResults`, optional 340 solar position parameters from a prior call to `solar_position`. 341 This can shorten runtime when calling multiple functions that use solar_position. 342 343 Returns 344 ------- 345 latitude : float 346 degrees of latitude 347 ''' 348 return solar_declination( *args, **kwargs )
Latitude of the subsolar point, degrees
Parameters
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- fast (bool):
see
solar_declination
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- latitude (float): degrees of latitude
350def solar_longitude( datetime, solar_pos=None ): 351 '''Longitude of the subsolar point, degrees 352 353 Parameters 354 ---------- 355 datetimeUTC : datetime-like or str 356 date and time. Include time zone or UTC will be assumed 357 solar_pos : `SolarPositionResults`, optional 358 solar position parameters from a prior call to `solar_position`. 359 This can shorten runtime when calling multiple functions that use solar_position. 360 361 Returns 362 ------- 363 longitude : float 364 degrees of longitude 365 ''' 366 # Convert to pandas Timestamp, if needed 367 datetimeUTC, tz_in = _to_datetime_utc(datetime) 368 369 # Longitude of subsolar point, degrees 370 # Equation of time will be added below 371 try: 372 # Treat as xarray.DataArray or pandas.Series 373 solar_lon = - 15 * ( datetimeUTC.dt.hour + 374 datetimeUTC.dt.minute / 60 + 375 datetimeUTC.dt.second / 3600 - 12 ) 376 except AttributeError: 377 solar_lon = - 15 * ( datetimeUTC.hour + 378 datetimeUTC.minute / 60 + 379 datetimeUTC.second / 3600 - 12 ) 380 381 # Add equation of time to the solar longitude, degrees 382 if solar_pos is None: 383 eot = equation_of_time( datetimeUTC, degrees=True ) 384 else: 385 eot = solar_pos.equation_of_time 386 solar_lon -= eot 387 388 return solar_lon
Longitude of the subsolar point, degrees
Parameters
- datetimeUTC (datetime-like or str): date and time. Include time zone or UTC will be assumed
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- longitude (float): degrees of longitude
390def solar_hour_angle( lon, datetime, solar_pos=None ): 391 '''Solar hour angle (degrees) for specified longitude, date and time 392 393 Hour angle is the angular displacement of the sun from the local meridian. 394 It is zero at local noon, negative in the morning, and positive is afternoon. 395 396 Parameters 397 ---------- 398 lon : float 399 longitude in degrees east 400 datetimeUTC : datetime-like or str 401 date and time. Include time zone or UTC will be assumed 402 solar_pos : `SolarPositionResults`, optional 403 solar position parameters from a prior call to `solar_position`. 404 This can shorten runtime when calling multiple functions that use solar_position. 405 406 Returns 407 ------- 408 ha : float 409 hour angle in degrees at the specified location and time 410 ''' 411 412 # Subsolar longitude, degrees 413 solar_lon = solar_longitude(datetime, solar_pos ) 414 415 # Hour angle, degrees 416 Ha = lon - solar_lon 417 418 return Ha
Solar hour angle (degrees) for specified longitude, date and time
Hour angle is the angular displacement of the sun from the local meridian. It is zero at local noon, negative in the morning, and positive is afternoon.
Parameters
- lon (float): longitude in degrees east
- datetimeUTC (datetime-like or str): date and time. Include time zone or UTC will be assumed
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- ha (float): hour angle in degrees at the specified location and time
420def equation_of_time( datetime, degrees=False, fast=False ): 421 '''Equation of time for specified date 422 423 Accounts for the solar day being slightly different from 24 hours 424 425 Parameters 426 ---------- 427 datetime : datetime-like or str 428 date and time. Include time zone or UTC will be assumed 429 degrees : bool (default=False) 430 If True, then return value in compass degrees 431 If False, then return value in minutes of an hour 432 fast : bool (default=False) 433 specifies whether to use a faster, but less accurate calculation 434 435 Returns 436 ------- 437 eot : float 438 equation of time on the specified date, degrees or minutes 439 ''' 440 # Convert to pandas Timestamp, if needed 441 datetime = _to_datetime(datetime) 442 443 # Determine whether to use the fast or accurate calculation 444 accurate = not fast 445 446 if accurate: 447 448 # Equation of time, minutes 449 junk, junk, eot, junk, junk = solar_position( datetime ) 450 451 else: 452 # Implements the "alternative equation" from Wikipedia, derived from 453 # https://web.archive.org/web/20120323231813/http://www.green-life-innovators.org/tiki-index.php?page=The%2BLatitude%2Band%2BLongitude%2Bof%2Bthe%2BSun%2Bby%2BDavid%2BWilliams 454 # When compared to the NREL SPA algorithm, differences reach are up to about 0.5 minute. 455 # Note: Leap years are not accounted for. 456 457 # Equation of time, accounts for the solar day differing slightly from 24 hr 458 try: 459 doy = datetime.dt.dayofyear 460 except AttributeError: 461 doy = datetime.dayofyear 462 W = 360 / 365.24 463 A = W * (doy+10) 464 B = A + 1.914 * np.sin( W * (doy-2) * pi180 ) 465 C = ( A - np.arctan2( np.tan(B*pi180), np.cos(23.44*pi180) ) / pi180 ) / 180 466 467 # Equation of time in minutes of an hour (1440 minutes per day) 468 eot = 720 * ( C - np.round(C) ) 469 470 # Equation of time, minutes -> degrees (360 degrees per day) 471 if degrees: 472 eot = eot / 60 * 360 / 24 473 474 return eot
Equation of time for specified date
Accounts for the solar day being slightly different from 24 hours
Parameters
- datetime (datetime-like or str): date and time. Include time zone or UTC will be assumed
- degrees (bool (default=False)): If True, then return value in compass degrees If False, then return value in minutes of an hour
- fast (bool (default=False)): specifies whether to use a faster, but less accurate calculation
Returns
- eot (float): equation of time on the specified date, degrees or minutes
476def solar_position( datetime ): 477 '''Compute position of sun (declination, right ascension, equation of time, distance)) on specified date 478 479 Calculations follow the NOAA solar calculator spreadsheet 480 Applicable to years 1900-2100. 481 482 Parameters 483 ---------- 484 date : datetime-like or str 485 date and time. Include time zone or UTC will be assumed 486 487 Returns 488 ------- 489 result : `SolarPositionResults` 490 declination, right ascension, equation of time, Earth-sun distance, datetime 491 ''' 492 # Ensure time is Timestamp in UTC 493 datetimeUTC, tz_in = _to_datetime_utc(datetime) 494 495 # Raise warning if any dates are outside date range 496 # recommended for orbital parameters used here 497 if np.logical_or( np.any( datetimeUTC < np.datetime64('1900-01-01') ), 498 np.any( datetimeUTC > np.datetime64('2100-01-01') ) ): 499 warnings.warn('Solar position accuracy declines for dates outside 1900-2100', \ 500 RuntimeWarning ) 501 502 # Number of days since 1 Jan 2000 503 NJD = datetimeUTC - np.datetime64('2000-01-01') 504 try: 505 NJD = NJD.dt.days \ 506 + NJD.dt.seconds / 86400. 507 except AttributeError: 508 NJD = NJD.days \ 509 + NJD.seconds / 86400. 510 511 # Julian day (since 12:00UTC 1 Jan 4713 BCE) 512 NJD += 2451544.50 513 514 # Julian century 515 JC = (NJD-2451545)/36525 516 517 # Earth orbital eccentricity, unitless 518 ec = 0.016708634 - JC*( 0.000042037 + 0.0000001267*JC ) 519 520 # Earth mean orbital obliquity, degree 521 mean_ob = 23 + ( 26 + ( (21.448 522 - JC * (46.815 523 + JC * (0.00059 - JC * 0.001813) ) ) )/60 )/60 524 525 # Earth true orbital obliquity, corrected for nutation, degree 526 ob = mean_ob + 0.00256 * np.cos( (125.04 - 1934.136*JC ) * pi180 ) 527 528 # Sun Mean ecliptic longitude, degree 529 mean_ec_lon = np.mod( 280.46646 + JC*( 36000.76983 + JC*0.0003032 ), 360 ) 530 531 # Sun Mean anomaly, degree 532 mean_anom = 357.52911 + JC*( 35999.05029 - 0.0001537*JC ) 533 534 # Sun Equation of center, degree 535 eq_center = np.sin(mean_anom*pi180) * (1.914602 - JC*( 0.004817 + 0.000014*JC )) \ 536 + np.sin(2*mean_anom*pi180) * (0.019993 - 0.000101*JC) \ 537 + np.sin(3*mean_anom*pi180) * 0.000289 538 539 # Sun True ecliptic longitude, degrees 540 true_ec_lon = mean_ec_lon + eq_center 541 542 # Sun True anomaly, degree 543 true_anom = mean_anom + eq_center 544 545 # Earth-Sun distance, AU 546 distance = (1.000001018 * (1-ec**2) ) / (1 + ec * np.cos( true_anom * pi180 )) 547 548 # Sun Apparent ecliptic longitude, corrected for nutation, degrees 549 ec_lon = true_ec_lon - 0.00569 - 0.00478 * np.sin( (125.04 - 1934.136*JC ) * pi180) 550 551 # Sun Right ascension, deg 552 right_ascension = np.arctan2( np.cos(ob*pi180) * np.sin(ec_lon*pi180), 553 np.cos(ec_lon*pi180) ) / pi180 554 555 # Sun Declination, deg 556 declination = np.arcsin( np.sin(ob*pi180) * np.sin(ec_lon*pi180) ) / pi180 557 558 # var y 559 vary = np.tan( ob/2 * pi180 )**2 560 561 # Equation of time, minutes 562 eot = vary * np.sin( 2 * mean_ec_lon * pi180) \ 563 - 2 * ec * np.sin( mean_anom * pi180 ) \ 564 + 4 * ec * vary * np.sin( mean_anom * pi180 ) * np.cos( 2 * mean_ec_lon * pi180) \ 565 - 0.5 * vary**2 * np.sin( 4 * mean_ec_lon * pi180) \ 566 - 1.25 * ec**2 * np.sin( 2 * mean_anom * pi180) 567 eot = eot * 4 / pi180 568 569 # Bundle results 570 result = SolarPositionResults(declination, 571 right_ascension, 572 eot, 573 distance, 574 datetimeUTC) 575 576 return result
Compute position of sun (declination, right ascension, equation of time, distance)) on specified date
Calculations follow the NOAA solar calculator spreadsheet Applicable to years 1900-2100.
Parameters
- date (datetime-like or str): date and time. Include time zone or UTC will be assumed
Returns
- result (
SolarPositionResults
): declination, right ascension, equation of time, Earth-sun distance, datetime
578def sun_times( lat, lon, datetime, tz_out=None, sza_sunrise=90.833, 579 fast=False, solar_pos=None ): 580 '''Compute times of sunrise, sunset, solar noon, and day length 581 582 Common options for solar zenith angle at sunrise 583 - 90.833 for first edge of sun rising, typical (0.567°) refraction (default) 584 - 90.267 for first edge of sun rising, no refraction 585 - 90 degrees for center of sun rising, no refraction 586 587 Note: horizon_zenith_angle can be used to compute a more accurate horizon location 588 for sites at elevation. 589 590 Parameters 591 ---------- 592 lat : float or ndarray 593 latitude in degrees 594 lon : float or ndarray 595 longitudes in degrees 596 datetime : datetime-like or str 597 datetime, provide a time zone or UTC will be assumed 598 tz_out : str, pytz.timezone, datetime.tzinfo, optional 599 timezone to be used for output times. 600 If None is provided, then result will be in same time zone as input or UTC 601 sza_sunrise : float (default=90.833) 602 Solar zenith angle at which sunrise and sunset are calculated, degrees 603 fast : bool (default=False) 604 Select a faster but less accurate calculation 605 solar_pos : `SolarPositionResults`, optional 606 solar position parameters from a prior call to `solar_position`. 607 This can shorten runtime when calling multiple functions that use solar_position. 608 609 Returns 610 ------- 611 result : `SunTimesResults` 612 times of sunrise, sunset, solar noon, and day length 613 ''' 614 # Convert to pandas Timestamp in UTC, if needed 615 datetimeUTC, tz_in = _to_datetime_utc(datetime) 616 617 # If no output timezone is specified, use the input time zone 618 if tz_out is None: 619 tz_out = tz_in 620 621 # Select fast or accurate calculation 622 accurate = not fast 623 624 # Solar declination (degrees) and equation of time (minutes) 625 if solar_pos is not None: 626 dec = solar_pos.declination 627 eot = solar_pos.eot 628 else: 629 if accurate: 630 dec, junk, eot, junk, junk = solar_position( datetimeUTC ) 631 else: 632 dec = solar_declination( datetimeUTC ) 633 eot = equation_of_time( datetimeUTC ) 634 635 # Sunrise hour angle, degree 636 # Degrees east of the local meridian where sun rises 637 ha_sunrise = np.arccos( np.cos(sza_sunrise*pi180) / 638 (np.cos(lat*pi180) * np.cos(dec*pi180)) 639 - np.tan(lat*pi180)*np.tan(dec*pi180) ) / pi180 640 641 # Solar noon, local standard time, day fraction 642 solar_noon = (720 - 4*lon - eot ) / 1440 643 644 # Sunrise and sunset, local standard time, day fraction 645 t_sunrise = solar_noon - 4 * ha_sunrise / 1440 646 t_sunset = solar_noon + 4 * ha_sunrise / 1440 647 648 # Midnight UTC 649 # datetimeUTC is in UTC but time-zone-naive 650 try: 651 # Series time objects 652 dateUTC = datetimeUTC.dt.floor('D') 653 except AttributeError: 654 # Scalar time objects 655 dateUTC = datetimeUTC.floor('D') 656 657 # Convert day fraction -> date time 658 solar_noon = dateUTC + solar_noon * pd.Timedelta( 1, 'day' ) 659 t_sunrise = dateUTC + t_sunrise * pd.Timedelta( 1, 'day' ) 660 t_sunset = dateUTC + t_sunset * pd.Timedelta( 1, 'day' ) 661 662 # Convert to output timezone, if any is provided 663 if tz_out is not None: 664 if isinstance(solar_noon,(xr.DataArray,np.ndarray)) or \ 665 isinstance(t_sunrise,(xr.DataArray,np.ndarray)): 666 # These types don't localize tz, but we can add offset to the tz-naive time 667 if hasattr(datetimeUTC,'tz_localize'): 668 # For scalar datetime, there is a single time offset, which we can add 669 utcoffset = np.timedelta64( datetimeUTC.tz_localize(tz_out).utcoffset() ) 670 solar_noon += utcoffset 671 t_sunrise += utcoffset 672 t_sunset += utcoffset 673 else: 674 # For Series datetime, there are potentially multiple offsets. We can only add one 675 unique_datetimeUTC = pd.DatetimeIndex(np.unique(datetimeUTC)) 676 unique_utcoffsets = np.unique( unique_datetimeUTC.tz_localize('UTC') \ 677 - unique_datetimeUTC.tz_localize(tz_out) ) 678 if len(unique_utcoffsets)==1: 679 utcoffset = unique_utcoffsets[0] 680 solar_noon += utcoffset 681 t_sunrise += utcoffset 682 t_sunset += utcoffset 683 else: 684 # We might be able to handle multiple offsets if we can 685 # determine which dimension of `solar_noon` and `t_sunrise` 686 # is the time dimension. 687 raise ValueError('Multiple timezone offsets not supported. ' 688 +'Request output in UTC or reduce number of input times.') 689 else: 690 try: 691 # Series time objects 692 solar_noon = solar_noon.dt.tz_localize('UTC')\ 693 .dt.tz_convert(tz_out) 694 t_sunrise = t_sunrise.dt.tz_localize('UTC')\ 695 .dt.tz_convert(tz_out) 696 t_sunset = t_sunset.dt.tz_localize('UTC')\ 697 .dt.tz_convert(tz_out) 698 except AttributeError: 699 # Scale time objects 700 solar_noon = solar_noon.tz_localize('UTC')\ 701 .tz_convert(tz_out) 702 t_sunrise = t_sunrise.tz_localize('UTC')\ 703 .tz_convert(tz_out) 704 t_sunset = t_sunset.tz_localize('UTC')\ 705 .tz_convert(tz_out) 706 707 # Sunlight duration, minutes 708 day_length = 8 * ha_sunrise * pd.Timedelta(1, 'minute') 709 710 result = SunTimesResults(t_sunrise, t_sunset, day_length, solar_noon, datetimeUTC) 711 712 return result
Compute times of sunrise, sunset, solar noon, and day length
Common options for solar zenith angle at sunrise
- 90.833 for first edge of sun rising, typical (0.567°) refraction (default)
- 90.267 for first edge of sun rising, no refraction
- 90 degrees for center of sun rising, no refraction
Note: horizon_zenith_angle can be used to compute a more accurate horizon location for sites at elevation.
Parameters
- lat (float or ndarray): latitude in degrees
- lon (float or ndarray): longitudes in degrees
- datetime (datetime-like or str): datetime, provide a time zone or UTC will be assumed
- tz_out (str, pytz.timezone, datetime.tzinfo, optional): timezone to be used for output times. If None is provided, then result will be in same time zone as input or UTC
- sza_sunrise (float (default=90.833)): Solar zenith angle at which sunrise and sunset are calculated, degrees
- fast (bool (default=False)): Select a faster but less accurate calculation
- solar_pos (
SolarPositionResults
, optional): solar position parameters from a prior call tosolar_position
. This can shorten runtime when calling multiple functions that use solar_position.
Returns
- result (
SunTimesResults
): times of sunrise, sunset, solar noon, and day length
714def horizon_zenith_angle( alt ): 715 '''Angle from the zenith to the horizon 716 717 The horizon is the locii of points where a line from the 718 observation location to the ellipsoid is tangent to the ellipsoid surface. 719 720 The altitude parameter should be the vertical distance 721 above the surrounding terrain that defines the horizon, 722 not necessarily the altitude above sea level or the altitude above ground level. 723 For example, on a mountain peak that is 4000 m above sea level and 724 1500 m above the surrounding plateau, the relevant altitude is 1500 m. 725 For an observer on the plateau, the relevant altitude is 0 m. 726 727 The implementation below assumes a spherical Earth. 728 Results using the WGS84 ellipsoid (see commented code below) 729 differ from the spherical case by << 1°. Terrain, 730 which is neglected here, has a larger effect on the horizon 731 location, so the simpler spherical calculation is appropriate. 732 733 Parameters 734 ---------- 735 lat : float or ndarray 736 latitude in degrees 737 alt : float or ndarray 738 altitude above surrounding terrain that defines the horizon, meters 739 740 Returns 741 ------- 742 hza : float or ndarray 743 horizon zenith angle in degrees 744 ''' 745 746 # WGS84 ellipsoid parameters 747 # semi-major radius, m 748 r_earth = 6378137.0 749 # ellipsoidal flattening, unitless 750 f = 1/298.257223563 751 752 # Horizon zenith angle, degrees (spherical earth) 753 hza = 180 - np.arcsin( r_earth / ( r_earth + alt ) ) / pi180 754 755 ## Ellipsoidal Earth 756 # # Eccentricity of ellipsoid 757 # ecc = f * (2-f) 758 # # Local (i.e. prime vertical) radius of curvature at latitude 759 # N = r_earth / np.sqrt( 1 - ecc**2 * np.sin(lat*pi180)**2 ) 760 # # Horizon zenith angle, degrees 761 # hza = 180 - np.arcsin( N / (N+alt) ) / pi180 762 763 return hza
Angle from the zenith to the horizon
The horizon is the locii of points where a line from the observation location to the ellipsoid is tangent to the ellipsoid surface.
The altitude parameter should be the vertical distance above the surrounding terrain that defines the horizon, not necessarily the altitude above sea level or the altitude above ground level. For example, on a mountain peak that is 4000 m above sea level and 1500 m above the surrounding plateau, the relevant altitude is 1500 m. For an observer on the plateau, the relevant altitude is 0 m.
The implementation below assumes a spherical Earth. Results using the WGS84 ellipsoid (see commented code below) differ from the spherical case by << 1°. Terrain, which is neglected here, has a larger effect on the horizon location, so the simpler spherical calculation is appropriate.
Parameters
- lat (float or ndarray): latitude in degrees
- alt (float or ndarray): altitude above surrounding terrain that defines the horizon, meters
Returns
- hza (float or ndarray): horizon zenith angle in degrees
765def refraction_angle( true_elevation_angle, pressure=101325., temperature_celsius=10. ): 766 '''Atmospheric refraction angle for light passing through Earth's atmosphere 767 768 The apparent locations in the sky of objects outsides Earth's atmosphere 769 differs from their true locations due to atmospheric refraction. 770 (e.g. sun and moon when they rise and set) 771 The apparent elevation of an object above the horizon is 772 apparent elevation angle = (true elevation angle) + (refraction angle) 773 774 The equations here are from Saemundsson/Bennett, whose calculations use 775 a typical vertical profile of atmospheric density (i.e. temperature and pressure). 776 The profiles can be rescaled to a particular surface temperature and pressure 777 to approximately account for varying atmospheric conditions. 778 Accurate refraction calculations should use fully specified vertical profile 779 of temperature and pressure, which cannot be done here. 780 781 Parameters 782 ---------- 783 true_elevation_angle : float 784 degrees above horizon of sun or other object 785 pressure : float (default=101325) 786 surface atmospheric pressure (Pa) 787 temperature_celsius : float (default=10) 788 surface atmospheric temperature (C) 789 790 Returns 791 ------- 792 angle : float 793 refraction angle in degrees. Value is zero when apparent elevation is below horizon 794 ''' 795 # Refraction angle, arcminutes 796 R = 1.02 / np.tan( ( true_elevation_angle + 10.3 / (true_elevation_angle + 5.11) ) * pi180 ) 797 # Account for temperature and pressure, arcminutes 798 R = R * pressure / 101325 * 283 / ( 273 + temperature_celsius ) 799 # Convert arcminutes -> degrees 800 R /= 60 801 802 # Result must be positive 803 R = np.maximum(R,0) 804 805 # Refraction defined only when the apparent elevation angle is positive 806 # Set refraction to zero when the apparent elevation is below horizon 807 refraction_angle = np.where( true_elevation_angle + R <= 0, 0, R) 808 809 return refraction_angle
Atmospheric refraction angle for light passing through Earth's atmosphere
The apparent locations in the sky of objects outsides Earth's atmosphere differs from their true locations due to atmospheric refraction. (e.g. sun and moon when they rise and set) The apparent elevation of an object above the horizon is apparent elevation angle = (true elevation angle) + (refraction angle)
The equations here are from Saemundsson/Bennett, whose calculations use a typical vertical profile of atmospheric density (i.e. temperature and pressure). The profiles can be rescaled to a particular surface temperature and pressure to approximately account for varying atmospheric conditions. Accurate refraction calculations should use fully specified vertical profile of temperature and pressure, which cannot be done here.
Parameters
- true_elevation_angle (float): degrees above horizon of sun or other object
- pressure (float (default=101325)): surface atmospheric pressure (Pa)
- temperature_celsius (float (default=10)): surface atmospheric temperature (C)
Returns
- angle (float): refraction angle in degrees. Value is zero when apparent elevation is below horizon
Namedtuple containing results of solar_position
All values are geocentric, strictly accurate for the center of the Earth, not a point on Earth's surface. The parallax angle from Earth's center to surface is 4e-5 degrees.
Attributes
- declination (float): position of the sun relative to Earth's equatorial plane, degrees
- right_ascension (float): position of the sun along Earth's equatorial plane relative to the vernal equinox, degrees
- equation_of_time (float): equation of time (minutes) between mean solar time and true solar time. Divide by 4 minutes per degree to obtain equation of time in degrees.
- distance (float): Earth-sun distance in AU (1 AU = 1.495978707e11 m)
- datetimeUTC (Timestamp, Series, or DataArray of numpy.datetime64): date and time input for calculations, UTC
Create new instance of SolarPositionResults(declination, right_ascension, equation_of_time, distance, datetimeUTC)
Inherited Members
- builtins.tuple
- index
- count
Namedtuple containing results of sun_times
Attributes
- sunrise (pandas.DatetimeIndex): sunrise time, UTC if not specified otherwise
- sunset (pandas.DatetimeIndex): sunset time, UTC if not specified otherwise
- day_length (pandas.Timedelta): duration of daylight
- solar_noon (pandas.DatetimeIndex): time of meridian transit, UTC if not specified otherwise
- datetimeUTC (Timestamp, Series, or DataArray of numpy.datetime64): date and time input for calculations, UTC
Create new instance of SunTimesResults(sunrise, sunset, day_length, solar_noon, datetimeUTC)
Inherited Members
- builtins.tuple
- index
- count
Constant $\pi/180$
159def solar_zenith_angle( lat, lon, datetime, 160 refraction=False, temperature=10., pressure=101325., 161 solar_pos=None ): 162 '''Solar zenith angle (degrees) for a given latitude, longitude, date and time. 163 164 Accounts for equation of time and (optionally) for atmospheric refraction. 165 Altitude of the observer is not accounted for, which can be important when the sun 166 is near the horizon. 167 168 Results are accurate to tenths of a degree, except where altitude is important 169 (< 20 degrees solar elevation) 170 171 Parameters 172 ---------- 173 lat : float or ndarray 174 latitude in degrees 175 lon : float or ndarray 176 longitudes in degrees 177 datetime : datetime-like or str 178 date and time. Include time zone or UTC will be assumed 179 refraction : bool (default=False) 180 specifies whether to account for atmospheric refraction 181 temperature : float or ndarray (default=10) 182 surface atmospheric temperature (Celsius), only used for refraction calculation 183 pressure : float or ndarray (default=101325) 184 surface atmospheric pressure (Pa), only used for refraction calculation 185 solar_pos : `SolarPositionResults`, optional 186 solar position parameters from a prior call to `solar_position`. 187 This can shorten runtime when calling multiple functions that use solar_position. 188 189 Returns 190 ------- 191 sza : float or ndarray 192 solar zenith angle in degrees at the designated locations and times 193 If refraction=False, this is the true solar zenith angle 194 If refraction=True, this is the apparent solar zenith angle 195 ''' 196 # Convert to pandas Timestamp, if needed 197 datetime = _to_datetime(datetime) 198 199 # Solar declination, degrees 200 if solar_pos is None: 201 dec = solar_declination( datetime ) 202 else: 203 dec = solar_pos.declination 204 205 # Hour angle, degrees 206 Ha = solar_hour_angle( lon, datetime, solar_pos ) 207 208 # True solar zenith angle, radians 209 sza = np.arccos( np.sin(lat*pi180) * np.sin(dec*pi180) + \ 210 np.cos(lat*pi180) * np.cos(dec*pi180) * np.cos(Ha*pi180) ) 211 212 # Convert radians -> degrees 213 sza /= pi180 214 215 if refraction: 216 # Subtract refraction angle (degrees) from zenith angle. 217 # SZA is always smaller due to refraction. 218 sza -= refraction_angle( 90-sza, pressure, temperature ) 219 220 return sza
Alias for solar_zenith_angle
63def solar_azimuth_angle( lat, lon, datetime, solar_pos=None ): 64 '''Solar azimuth angle (degrees) for a latitude, longitude, date and time 65 66 SAA is degrees clockwise from north. 67 68 Parameters 69 ---------- 70 lat : float or ndarray 71 latitude in degrees 72 lon : float or ndarray 73 longitudes in degrees 74 datetime : datetime-like or str 75 date and time. Include time zone or UTC will be assumed 76 solar_pos : `SolarPositionResults`, optional 77 solar position parameters from a prior call to `solar_position`. 78 This can shorten runtime when calling multiple functions that use solar_position. 79 80 Returns 81 ------- 82 saa : float or ndarray 83 solar azimuth angle in degrees (clockwise from north) 84 ''' 85 # Convert to pandas Timestamp, if needed 86 datetime = _to_datetime(datetime) 87 88 # Subsolar point, latitude longitude, degrees 89 solar_lat = solar_latitude( datetime, solar_pos=solar_pos ) 90 solar_lon = solar_longitude( datetime, solar_pos=solar_pos ) 91 92 # Vector pointing toward sun 93 x = np.cos( solar_lat * pi180 ) * np.sin( (solar_lon - lon) * pi180 ) 94 y = np.cos( lat*pi180 ) * np.sin( solar_lat*pi180 ) \ 95 - np.sin( lat*pi180 ) * np.cos( solar_lat*pi180 ) \ 96 * np.cos( (solar_lon - lon) * pi180 ) 97 98 # Azimuth angle from north, degrees 99 saa = np.arctan2( x, y ) / pi180 100 101 # Change range [-180,180] to [0,360] 102 return np.mod( saa+360, 360 )
Alias for solar_azimuth_angle
104def solar_elevation_angle( lat, lon, datetime, alt=0, 105 refraction=False, temperature=10., pressure=101325., 106 solar_pos=None ): 107 '''Solar elevation angle (degrees) above the horizon 108 109 The altitude parameter should be the vertical distance 110 above the surrounding terrain that defines the horizon, 111 not necessarily the altitude above sea level or the altitude above ground level. 112 For example, on a mountain peak that is 4000 m above sea level and 113 1500 m above the surrounding plateau, the relevant altitude is 1500 m. 114 For an observer on the plateau, the relevant altitude is 0 m. 115 116 See documentation for `solar_zenith_angle` and `horizon_zenith_angle`. 117 118 Parameters 119 ---------- 120 lat : float or ndarray 121 latitude in degrees 122 lon : float or ndarray 123 longitudes in degrees 124 datetime : datetime-like or str 125 date and time. Include time zone or UTC will be assumed 126 alt : float or ndarray (default=0) 127 altitude above surrounding terrain that defines the horizon, meters 128 refraction : bool (default=False) 129 specifies whether to account for atmospheric refraction 130 temperature : float or ndarray (default=10) 131 surface atmospheric temperature (Celsius), only used for refraction calculation 132 pressure : float or ndarray (default=101325) 133 surface atmospheric pressure (Pa), only used for refraction calculation 134 solar_pos : `SolarPositionResults`, optional 135 solar position parameters from a prior call to `solar_position`. 136 This can shorten runtime when calling multiple functions that use solar_position. 137 138 Returns 139 ------- 140 sea : float or ndarray 141 solar elevation angle in degrees at the designated locations and times 142 If refraction=False, this is the true solar elevation angle 143 If refraction=True, this is the apparent solar elevation angle 144 145 ''' 146 147 if refraction and np.any(alt): 148 warnings.warn( 'Atmospheric refraction is calculated for surface conditions, ' 149 + 'but an altitude above the surface was specified', 150 category=UserWarning, 151 stacklevel=2 ) 152 153 sea = horizon_zenith_angle( alt ) \ 154 - solar_zenith_angle( lat, lon, datetime, refraction, temperature, pressure, 155 solar_pos=solar_pos ) 156 157 return sea
Alias for solar_elevation_angle
420def equation_of_time( datetime, degrees=False, fast=False ): 421 '''Equation of time for specified date 422 423 Accounts for the solar day being slightly different from 24 hours 424 425 Parameters 426 ---------- 427 datetime : datetime-like or str 428 date and time. Include time zone or UTC will be assumed 429 degrees : bool (default=False) 430 If True, then return value in compass degrees 431 If False, then return value in minutes of an hour 432 fast : bool (default=False) 433 specifies whether to use a faster, but less accurate calculation 434 435 Returns 436 ------- 437 eot : float 438 equation of time on the specified date, degrees or minutes 439 ''' 440 # Convert to pandas Timestamp, if needed 441 datetime = _to_datetime(datetime) 442 443 # Determine whether to use the fast or accurate calculation 444 accurate = not fast 445 446 if accurate: 447 448 # Equation of time, minutes 449 junk, junk, eot, junk, junk = solar_position( datetime ) 450 451 else: 452 # Implements the "alternative equation" from Wikipedia, derived from 453 # https://web.archive.org/web/20120323231813/http://www.green-life-innovators.org/tiki-index.php?page=The%2BLatitude%2Band%2BLongitude%2Bof%2Bthe%2BSun%2Bby%2BDavid%2BWilliams 454 # When compared to the NREL SPA algorithm, differences reach are up to about 0.5 minute. 455 # Note: Leap years are not accounted for. 456 457 # Equation of time, accounts for the solar day differing slightly from 24 hr 458 try: 459 doy = datetime.dt.dayofyear 460 except AttributeError: 461 doy = datetime.dayofyear 462 W = 360 / 365.24 463 A = W * (doy+10) 464 B = A + 1.914 * np.sin( W * (doy-2) * pi180 ) 465 C = ( A - np.arctan2( np.tan(B*pi180), np.cos(23.44*pi180) ) / pi180 ) / 180 466 467 # Equation of time in minutes of an hour (1440 minutes per day) 468 eot = 720 * ( C - np.round(C) ) 469 470 # Equation of time, minutes -> degrees (360 degrees per day) 471 if degrees: 472 eot = eot / 60 * 360 / 24 473 474 return eot
Alias for equation_of_time
271def solar_declination( datetime, fast=False, solar_pos=None ): 272 '''Calculate solar declination (degrees) for specified date 273 274 Parameters 275 ---------- 276 datetime : datetime-like or str 277 date and time. Include time zone or UTC will be assumed 278 fast : bool (default=False) 279 Specifies using a faster but less accurate calculation 280 solar_pos : `SolarPositionResults`, optional 281 solar position parameters from a prior call to `solar_position`. 282 This can shorten runtime when calling multiple functions that use solar_position. 283 284 Returns 285 ------- 286 dec : float 287 solar declination in degrees at the specified date 288 ''' 289 # Convert to pandas Timestamp, if needed 290 datetime = _to_datetime(datetime) 291 292 # Select the accurate or fast calculation 293 accurate = not fast 294 295 if solar_pos is not None: 296 dec = solar_pos.declination 297 298 else: 299 if accurate: 300 301 # Solar declination, degrees 302 dec, junk, junk, junk, junk = solar_position( datetime ) 303 304 else: 305 # The fast method implements 306 # Eq. 9.68-9.72 from M.Z. Jacobson, Fundamentals of Atmospheric Modeling 307 308 # Number of days since beginning of 2000 309 NJD = datetime - np.datetime64('2000-01-01') 310 try: 311 NJD = NJD.dt.days 312 except AttributeError: 313 NJD = NJD.days 314 315 # Obliquity, degrees 316 ob = 23.439 - 4e-7 * NJD 317 318 # Parameters for ecliptic, degrees 319 gm = 357.528 + 0.9856003 * NJD 320 lm = 280.460 + 0.9856474 * NJD 321 322 # Ecliptic longitude of sun, degrees 323 ec = lm + 1.915 * np.sin( gm * pi180 ) + 0.020 * np.sin( 2 * gm * pi180 ) 324 325 #Solar declination, degrees 326 dec = np.arcsin( np.sin( ob * pi180 ) * np.sin( ec * pi180 ) ) / pi180 327 328 return dec
Alias for solar_declination