Package spharm :: Module spharm

Source Code for Module spharm.spharm

   1  """ 
   2  Introduction 
   3  ============ 
   4       
   5  This module provides a python interface to the NCAR  
   6  U{SPHEREPACK<http://www.scd.ucar.edu/softlib/SPHERE.html>} library.   
   7  It is not a one-to-one wrapper for the SPHEREPACK routines, rather 
   8  it provides a simple interface oriented toward working with 
   9  atmospheric general circulation model (GCM) data. 
  10   
  11  Requirements 
  12  ============ 
  13   - U{numpy<http://numeric.scipy.org>}, and a fortran compiler 
  14   supported by numpy.f2py. 
  15   
  16  Installation 
  17  ============ 
  18   - U{Download<http://code.google.com/p/pyspharm/downloads/list>} module source, 
  19   untar. 
  20   - run C{python setup.py install} (as root if necessary). 
  21   The SPHERPACK fortran source files will be downloaded automatically by the 
  22   setup.py script, since the SPHEREPACK license prohibits redistribution. 
  23   To specify the fortran compiler to use (e.g. g95) run 
  24   C{python setup.py config_fc --fcompiler=g95 install}. C{f2py -c --help-fcompiler} 
  25   will show you what fortran compilers are available. 
  26   
  27  Usage 
  28  ===== 
  29   
  30  >>> import spharm 
  31  >>> x=spharm.Spharmt(144,72,rsphere=8e6,gridtype='gaussian',legfunc='computed') 
  32   
  33  creates a class instance for spherical harmonic calculations on a 144x72 
  34  gaussian grid on a sphere with radius 8000 km. The associated legendre  
  35  functions are recomputed on the fly (instead of pre-computed and stored). 
  36  Default values of rsphere, gridtype and legfunc are 6.3712e6, 'regular' 
  37  and 'stored'. Real-world examples are included in the source distribution. 
  38   
  39  Class methods 
  40  ============= 
  41   - grdtospec: grid to spectral transform (spherical harmonic analysis). 
  42   - spectogrd: spectral to grid transform (spherical harmonic synthesis). 
  43   - getuv:  compute u and v winds from spectral coefficients of vorticity 
  44   and divergence. 
  45   - getvrtdivspec: get spectral coefficients of vorticity and divergence 
  46   from u and v winds. 
  47   - getgrad: compute the vector gradient given spectral coefficients. 
  48   - getpsichi: compute streamfunction and velocity potential from winds. 
  49   - specsmooth:  isotropic spectral smoothing. 
  50   
  51  Functions 
  52  ========= 
  53   - regrid:  spectral re-gridding, with optional spectral smoothing and/or 
  54   truncation. 
  55   - gaussian_lats_wts: compute gaussian latitudes and weights. 
  56   - getspecindx: compute indices of zonal wavenumber and degree 
  57   for complex spherical harmonic coefficients. 
  58   - legendre: compute associated legendre functions. 
  59   - getgeodesicpts: computes the points on the surface of the sphere 
  60   corresponding to a twenty-sided (icosahedral) geodesic. 
  61   - specintrp: spectral interpolation to an arbitrary point on the sphere. 
  62   
  63  Conventions 
  64  =========== 
  65       
  66  The gridded data is assumed to be oriented such that i=1 is the  
  67  Greenwich meridian and j=1 is the northernmost point. Grid indices 
  68  increase eastward and southward. If nlat is odd the equator is included. 
  69  If nlat is even the equator will lie half way between points nlat/2 
  70  and (nlat/2)+1. nlat must be at least 3. For regular grids  
  71  (gridtype='regular') the poles will be included when nlat is odd. 
  72  The grid increment in longitude is 2*pi/nlon radians. For example, 
  73  nlon = 72 for a five degree grid. nlon must be greater than or  
  74  equal to 4. The efficiency of the computation is improved when nlon 
  75  is a product of small prime numbers. 
  76   
  77  The spectral data is assumed to be in a complex array of dimension 
  78  (ntrunc+1)*(ntrunc+2)/2. ntrunc is the triangular truncation limit 
  79  (ntrunc = 42 for T42). ntrunc must be <= nlat-1. Coefficients are 
  80  ordered so that first (nm=0) is m=0,n=0, second is m=0,n=1, 
  81  nm=ntrunc is m=0,n=ntrunc, nm=ntrunc+1 is m=1,n=1, etc.  
  82  The values of m (degree) and n (order) as a function of the index  
  83  nm are given by the arrays indxm, indxn returned by getspecindx. 
  84   
  85  The associated legendre polynomials are normalized so that the 
  86  integral (pbar(n,m,theta)**2)*sin(theta) on the interval theta=0 to pi 
  87  is 1, where pbar(m,n,theta)=sqrt((2*n+1)*factorial(n-m)/(2*factorial(n+m)))* 
  88  sin(theta)**m/(2**n*factorial(n)) times the (n+m)th derivative of 
  89  (x**2-1)**n with respect to x=cos(theta). 
  90  theta = pi/2 - phi, where phi is latitude and theta is colatitude. 
  91  Therefore, cos(theta) = sin(phi) and sin(theta) = cos(phi). 
  92  Note that pbar(0,0,theta)=sqrt(2)/2, and pbar(1,0,theta)=.5*sqrt(6)*sin(lat). 
  93   
  94  The default grid type is regular (equally spaced latitude points).   
  95  Set gridtype='gaussian' when creating a class instance 
  96  for gaussian latitude points. 
  97   
  98  Quantities needed to compute spherical harmonics are precomputed and stored 
  99  when the class instance is created with legfunc='stored' (the default). 
 100  If legfunc='computed', they are recomputed on the fly on each method call. 
 101  The storage requirements for legfunc="stored" increase like nlat**2, while 
 102  those for legfunc='stored' increase like nlat**3.  However, for 
 103  repeated method invocations on a single class instance, legfunc="stored" 
 104  will always be faster. 
 105   
 106  @contact: U{Jeff Whitaker<mailto:jeffrey.s.whitaker@noaa.gov>} 
 107   
 108  @version: 1.0.5       
 109   
 110  @license: Permission to use, copy, modify, and distribute this software and its 
 111  documentation for any purpose and without fee is hereby granted, 
 112  provided that the above copyright notice appear in all copies and that 
 113  both that copyright notice and this permission notice appear in 
 114  supporting documentation. 
 115  THE AUTHOR DISCLAIMS ALL WARRANTIES WITH REGARD TO THIS SOFTWARE, 
 116  INCLUDING ALL IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS. IN NO 
 117  EVENT SHALL THE AUTHOR BE LIABLE FOR ANY SPECIAL, INDIRECT OR 
 118  CONSEQUENTIAL DAMAGES OR ANY DAMAGES WHATSOEVER RESULTING FROM LOSS OF 
 119  USE, DATA OR PROFITS, WHETHER IN AN ACTION OF CONTRACT, NEGLIGENCE OR 
 120  OTHER TORTIOUS ACTION, ARISING OUT OF OR IN CONNECTION WITH THE USE OR 
 121  PERFORMANCE OF THIS SOFTWARE. 
 122  """ 
 123  import _spherepack, numpy, math, sys 
 124   
 125  # define a list of instance variables that cannot be rebound 
 126  # or unbound. 
 127  _private_vars = ['nlon','nlat','gridtype','legfunc','rsphere'] 
 128  __version__ = '1.0.5' 
 129   
130 -class Spharmt:
131 """ 132 spherical harmonic transform class. 133 134 @ivar nlat: number of latitudes (set when class instance is created, 135 cannot be changed). 136 137 @ivar nlon: number of longitudes (set when class instance is created, 138 cannot be changed). 139 140 @ivar rsphere: The radius of the sphere in meters (set when class 141 instance is created, cannot be changed). 142 143 @ivar legfunc: 'stored' or 'computed'. If 'stored', 144 associated legendre functions are precomputed and stored when the 145 class instance is created. If 'computed', associated 146 legendre functions are computed on the fly when transforms are 147 requested. Set when class instance is created, cannot be changed. 148 149 @ivar gridtype: 'regular' (equally spaced in longitude and latitude) 150 or 'gaussian' (equally spaced in longitude, latitudes located at 151 roots of ordinary Legendre polynomial of degree nlat). Set when class 152 instance is created, cannot be changed. 153 """ 154
155 - def __setattr__(self, key, val):
156 """ 157 prevent modification of read-only instance variables. 158 """ 159 if self.__dict__.has_key(key) and key in _private_vars: 160 raise AttributeError, 'Attempt to rebind read-only instance variable '+key 161 else: 162 self.__dict__[key] = val
163
164 - def __delattr__(self, key):
165 """ 166 prevent deletion of read-only instance variables. 167 """ 168 if self.__dict__.has_key(key) and key in _private_vars: 169 raise AttributeError, 'Attempt to unbind read-only instance variable '+key 170 else: 171 del self.__dict__[key]
172
173 - def __init__(self, nlon, nlat, rsphere=6.3712e6, gridtype='regular', legfunc='stored'):
174 """ 175 create a Spharmt class instance. 176 177 @param nlon: Number of longitudes. The grid must be oriented from 178 east to west, with the first point at the Greenwich meridian 179 and the last point at 360-delta degrees east 180 (where delta = 360/nlon degrees). Must be >= 4. Transforms will 181 be faster when nlon is the product of small primes. 182 183 @param nlat: Number of latitudes. The grid must be oriented from north 184 to south. If nlat is odd the equator is included. 185 If nlat is even the equator will lie half way between points 186 points nlat/2 and (nlat/2)+1. Must be >=3. 187 188 @keyword rsphere: The radius of the sphere in meters. 189 Default 6371200 (the value for Earth). 190 191 @keyword legfunc: 'stored' (default) or 'computed'. If 'stored', 192 associated legendre functions are precomputed and stored when the 193 class instance is created. This uses O(nlat**3) memory, but 194 speeds up the spectral transforms. If 'computed', associated 195 legendre functions are computed on the fly when transforms are 196 requested. This uses O(nlat**2) memory, but slows down the spectral 197 transforms a bit. 198 199 @keyword gridtype: 'regular' (default) or 'gaussian'. Regular grids 200 will include the poles and equator if nlat is odd. Gaussian 201 grids never include the poles, but will include the equator if 202 nlat is odd. 203 204 """ 205 # sanity checks. 206 if rsphere > 0.0: 207 self.rsphere= rsphere 208 else: 209 msg = 'Spharmt.__init__ illegal value of rsphere (%s) - must be postitive' % (rsphere) 210 raise ValueError, msg 211 if nlon > 3: 212 self.nlon = nlon 213 else: 214 msg = 'Spharmt.__init__ illegal value of nlon (%s) - must be at least 4' % (nlon,) 215 raise ValueError, msg 216 if nlat > 2: 217 self.nlat = nlat 218 else: 219 msg = 'Spharmt.__init__ illegal value of nlat (%s) - must be at least 3' % (nlat,) 220 raise ValueError, msg 221 if gridtype != 'regular' and gridtype != 'gaussian': 222 msg = 'Spharmt.__init__ illegal value of gridtype (%s) - must be either "gaussian" or "regular"' % gridtype 223 raise ValueError, msg 224 else: 225 self.gridtype = gridtype 226 227 if legfunc != 'computed' and legfunc != 'stored': 228 msg = 'Spharmt.__init__ illegal value of legfunc (%s) - must be either "computed" or "stored"' % legfunc 229 raise ValueError, msg 230 else: 231 self.legfunc = legfunc 232 233 if nlon%2: # nlon is odd 234 n1 = min(nlat, (nlon + 1)/2) 235 else: 236 n1 = min(nlat, (nlon + 2)/2) 237 if nlat%2: # nlat is odd 238 n2 = (nlat + 1)/2 239 else: 240 n2 = nlat/2 241 242 if gridtype == 'regular': 243 if legfunc == 'stored': 244 lshaes = (n1*n2*(nlat + nlat - n1 + 1))/2 + nlon + 15 245 lwork = 5*nlat*n2 + 3*((n1 - 2)*(nlat + nlat - n1 -1))/2 246 wshaes, ierror = _spherepack.shaesi(nlat, nlon, lshaes, lwork, nlat+1) 247 if ierror != 0: 248 msg = 'In return from call to shaesi in Spharmt.__init__ ierror = %d' % ierror 249 raise ValueError, msg 250 self.wshaes = wshaes 251 lshses = lshaes 252 wshses, ierror = _spherepack.shsesi(nlat, nlon, lshses, lwork, nlat+1) 253 if ierror != 0: 254 msg = 'In return from call to shsesi in Spharmt.__init__ ierror = %d' % ierror 255 raise ValueError, msg 256 self.wshses = wshses 257 lvhaes = n1*n2*(nlat + nlat - n1 + 1) + nlon + 15 258 lwork = 3*(max(n1 -2,0)*(nlat + nlat - n1 - 1))/2 + 5*n2*nlat 259 wvhaes, ierror = _spherepack.vhaesi(nlat, nlon, lvhaes, lwork, 2*(nlat+1)) 260 if ierror != 0: 261 msg = 'In return from call to vhaesi in Spharmt.__init__ ierror = %d' % ierror 262 raise ValueError, msg 263 self.wvhaes = wvhaes 264 lwork = 3*(max(n1 - 2,0)*(nlat + nlat - n1 -1))/2 + 5*n2*nlat 265 lvhses = n1*n2*(nlat + nlat - n1 + 1) + nlon + 15 266 wvhses, ierror = _spherepack.vhsesi(nlat,nlon,lvhses,lwork,2*(nlat+1)) 267 if ierror != 0: 268 msg = 'In return from call to vhsesi in Spharmt.__init__ ierror = %d' % ierror 269 raise ValueError, msg 270 self.wvhses = wvhses 271 else: 272 lshaec = 2*nlat*n2+3*((n1-2)*(nlat+nlat-n1-1))/2+nlon+15 273 wshaec, ierror = _spherepack.shaeci(nlat, nlon, lshaec, 2*(nlat+1)) 274 if ierror != 0: 275 msg = 'In return from call to shaeci in Spharmt.__init__ ierror = %d' % ierror 276 raise ValueError, msg 277 self.wshaec = wshaec 278 lshsec = lshaec 279 wshsec, ierror = _spherepack.shseci(nlat, nlon, lshsec, 2*(nlat+1)) 280 if ierror != 0: 281 msg = 'In return from call to shseci in Spharmt.__init__ ierror = %d' % ierror 282 raise ValueError, msg 283 self.wshsec = wshsec 284 lvhaec = 4*nlat*n2+3*max(n1-2,0)*(2*nlat-n1-1)+nlon+15 285 wvhaec, ierror = _spherepack.vhaeci(nlat, nlon, lvhaec, 2*(nlat+1)) 286 if ierror != 0: 287 msg = 'In return from call to vhaeci in Spharmt.__init__ ierror = %d' % ierror 288 raise ValueError, msg 289 self.wvhaec = wvhaec 290 lvhsec = lvhaec 291 wvhsec, ierror = _spherepack.vhseci(nlat, nlon, lvhsec, 2*(nlat+1)) 292 if ierror != 0: 293 msg = 'In return from call to vhseci in Spharmt.__init__ ierror = %d' % ierror 294 raise ValueError, msg 295 self.wvhsec = wvhsec 296 297 298 elif gridtype == 'gaussian': 299 if legfunc == 'stored': 300 lshags = nlat*(3*(n1 + n2) - 2) + (n1 - 1)*(n2*(2*nlat - n1) - 3*n1)/2 + nlon + 15 301 lwork = 4*nlat*(nlat + 2) + 2 302 ldwork = nlat*(nlat + 4) 303 wshags, ierror = _spherepack.shagsi(nlat, nlon, lshags, lwork, ldwork) 304 if ierror != 0: 305 msg = 'In return from call to shagsi in Spharmt.__init__ ierror = %d' % ierror 306 raise ValueError, msg 307 self.wshags = wshags 308 lshsgs = lshags 309 wshsgs, ierror = _spherepack.shsgsi(nlat, nlon, lshsgs, lwork, ldwork) 310 if ierror != 0: 311 msg = 'In return from call to shsgsi in Spharmt.__init__ ierror = %d' % ierror 312 raise ValueError, msg 313 self.wshsgs = wshsgs 314 lvhags = (nlat +1)*(nlat + 1)*nlat/2 +nlon + 15 315 ldwork = (3*nlat*(nlat + 3) + 2)/2 316 wvhags, ierror = _spherepack.vhagsi(nlat, nlon, lvhags, ldwork) 317 if ierror != 0: 318 msg = 'In return from call to vhagsi in Spharmt.__init__ ierror = %d' % ierror 319 raise ValueError, msg 320 self.wvhags = wvhags 321 lvhsgs = n1*n2*(nlat + nlat - n1 +1) + nlon + 15 + 2*nlat 322 ldwork = (3*nlat*(nlat + 3) + 2)/2 323 wvhsgs, ierror = _spherepack.vhsgsi(nlat, nlon, lvhsgs, ldwork) 324 if ierror != 0: 325 msg = 'In return from call to vhsgsi in Spharmt.__init__ ierror = %d' % ierror 326 raise ValueError, msg 327 self.wvhsgs = wvhsgs 328 else: 329 lshagc = nlat*(2*n2+3*n1-2)+3*n1*(1-n1)/2+nlon+15 330 wshagc, ierror = _spherepack.shagci(nlat, nlon, lshagc, nlat*(nlat+4)) 331 if ierror != 0: 332 msg = 'In return from call to shagci in Spharmt.__init__ ierror = %d' % ierror 333 raise ValueError, msg 334 self.wshagc = wshagc 335 lshsgc = lshagc 336 wshsgc, ierror = _spherepack.shsgci(nlat, nlon, lshsgc, nlat*(nlat+4)) 337 if ierror != 0: 338 msg = 'In return from call to shsgci in Spharmt.__init__ ierror = %d' % ierror 339 raise ValueError, msg 340 self.wshsgc = wshsgc 341 lvhagc = 4*nlat*n2+3*max(n1-2,0)*(2*nlat-n1-1)+nlon+n2+15 342 ldwork = 2*nlat*(nlat+1)+1 343 wvhagc, ierror = _spherepack.vhagci(nlat, nlon, lvhagc, ldwork) 344 if ierror != 0: 345 msg = 'In return from call to vhagci in Spharmt.__init__ ierror = %d' % ierror 346 raise ValueError, msg 347 self.wvhagc = wvhagc 348 lvhsgc = 4*nlat*n2+3*max(n1-2,0)*(2*nlat-n1-1)+nlon+15 349 wvhsgc, ierror = _spherepack.vhsgci(nlat, nlon, lvhsgc, ldwork) 350 if ierror != 0: 351 msg = 'In return from call to vhsgci in Spharmt.__init__ ierror = %d' % ierror 352 raise ValueError, msg 353 self.wvhsgc = wvhsgc
354
355 - def grdtospec(self, datagrid, ntrunc=None):
356 """ 357 grid to spectral transform (spherical harmonic analysis). 358 359 @param datagrid: rank 2 or 3 numpy float32 array with shape (nlat,nlon) or 360 (nlat,nlon,nt), where nt is the number of grids to be transformed. If 361 datagrid is rank 2, nt is assumed to be 1. 362 363 @keyword ntrunc: optional spectral truncation limit. 364 (default self.nlat-1) 365 366 @return: C{B{dataspec}} - rank 1 or 2 numpy complex array with shape 367 (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing 368 complex spherical harmonic coefficients resulting from the spherical 369 harmonic analysis of datagrid. 370 """ 371 372 # check that datagrid is rank 2 or 3 with size (self.nlat, self.nlon) or 373 # (self.nlat, self.nlon, nt) where nt is number of grids to transform. 374 375 if len(datagrid.shape) > 3: 376 msg = 'grdtospec needs a rank two or three array, got %d' % (len(datagrid.shape),) 377 raise ValueError, msg 378 379 if datagrid.shape[0] != self.nlat or datagrid.shape[1] != self.nlon: 380 msg = 'grdtospec needs an array of size %d by %d, got %d by %d' % (self.nlat, self.nlon, datagrid.shape[0], datagrid.shape[1],) 381 raise ValueError, msg 382 383 # check ntrunc. 384 385 if ntrunc is None: 386 ntrunc = self.nlat-1 387 388 if ntrunc < 0 or ntrunc+1 > datagrid.shape[0]: 389 msg = 'ntrunc must be between 0 and %d' % (datagrid.shape[0]-1,) 390 raise ValueError, msg 391 392 393 nlat = self.nlat 394 nlon = self.nlon 395 if nlat%2: # nlat is odd 396 n2 = (nlat + 1)/2 397 else: 398 n2 = nlat/2 399 400 if len(datagrid.shape) == 2: 401 nt = 1 402 datagrid = numpy.reshape(datagrid, (nlat,nlon,1)) 403 else: 404 nt = datagrid.shape[2] 405 406 if self.gridtype == 'regular': 407 408 # do grid to spectral transform. 409 410 if self.legfunc == 'stored': 411 lwork = (nt+1)*nlat*nlon 412 a,b,ierror = _spherepack.shaes(datagrid,self.wshaes,lwork) 413 if ierror != 0: 414 msg = 'In return from call to shaes in Spharmt.grdtospec ierror = %d' % ierror 415 raise ValueError, msg 416 else: 417 lwork = nlat*(nt*nlon+max(3*n2,nlon)) 418 a,b,ierror = _spherepack.shaec(datagrid,self.wshaec,lwork) 419 if ierror != 0: 420 msg = 'In return from call to shaec in Spharmt.grdtospec ierror = %d' % ierror 421 422 raise ValueError, msg 423 424 # gaussian grid. 425 426 elif self.gridtype == 'gaussian': 427 428 # do grid to spectral transform. 429 430 if self.legfunc == 'stored': 431 lwork = nlat*nlon*(nt+1) 432 a,b,ierror = _spherepack.shags(datagrid,self.wshags,lwork) 433 if ierror != 0: 434 msg = 'In return from call to shags in Spharmt.grdtospec ierror = %d' % ierror 435 raise ValueError, msg 436 else: 437 lwork = nlat*(nlon*nt+max(3*n2,nlon)) 438 a,b,ierror = _spherepack.shagc(datagrid,self.wshagc,lwork) 439 if ierror != 0: 440 msg = 'In return from call to shagc in Spharmt.grdtospec ierror = %d' % ierror 441 raise ValueError, msg 442 443 # convert 2d real and imag spectral arrays into 1d complex array. 444 445 dataspec = _spherepack.twodtooned(a,b,ntrunc) 446 447 if nt==1: 448 return numpy.squeeze(dataspec) 449 else: 450 return dataspec
451
452 - def spectogrd(self, dataspec):
453 454 """ 455 spectral to grid transform (spherical harmonic synthesis). 456 457 @param dataspec: rank 1 or 2 numpy complex array with shape 458 (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing 459 complex spherical harmonic coefficients (where ntrunc is the 460 triangular truncation limit and nt is the number of spectral arrays 461 to be transformed). If dataspec is rank 1, nt is assumed to be 1. 462 463 @return: C{B{datagrid}} - rank 2 or 3 numpy float32 array with shape 464 (nlat,nlon) or (nlat,nlon,nt) containing the gridded data resulting from 465 the spherical harmonic synthesis of dataspec. 466 """ 467 468 # make sure dataspec is rank 1 or 2. 469 470 if len(dataspec.shape) > 2: 471 msg = 'spectogrd needs a rank one or two array, got %d' % (len(dataspec.shape),) 472 raise ValueError, msg 473 474 nlat = self.nlat 475 nlon = self.nlon 476 if nlat%2: # nlat is odd 477 n2 = (nlat + 1)/2 478 else: 479 n2 = nlat/2 480 481 if len(dataspec.shape) == 1: 482 nt = 1 483 dataspec = numpy.reshape(dataspec, (dataspec.shape[0],1)) 484 else: 485 nt = dataspec.shape[1] 486 487 ntrunc = int(-1.5 + 0.5*math.sqrt(9.-8.*(1.-dataspec.shape[0]))) 488 if ntrunc > nlat-1: 489 msg = 'ntrunc too large - can be max of %d, got %d' % (nlat-1,ntrunc) 490 raise ValueError, msg 491 492 a, b = _spherepack.onedtotwod(dataspec,nlat) 493 494 # regular grid. 495 496 if self.gridtype == 'regular': 497 498 # do spectral to grid transform. 499 500 if self.legfunc == 'stored': 501 lwork = (nt+1)*nlat*nlon 502 datagrid, ierror = _spherepack.shses(nlon,a,b,self.wshses,lwork) 503 if ierror != 0: 504 msg = 'In return from call to shses in Spharmt.spectogrd ierror = %d' % ierror 505 raise ValueError, msg 506 else: 507 lwork = nlat*(nt*nlon+max(3*n2,nlon)) 508 datagrid, ierror = _spherepack.shsec(nlon,a,b,self.wshsec,lwork) 509 if ierror != 0: 510 msg = 'In return from call to shsec in Spharmt.spectogrd ierror = %d' % ierror 511 raise ValueError, msg 512 513 # gaussian grid. 514 515 elif self.gridtype == 'gaussian': 516 517 # do spectral to grid transform. 518 519 if self.legfunc == 'stored': 520 lwork = nlat*nlon*(nt+1) 521 datagrid, ierror = _spherepack.shsgs(nlon,a,b,self.wshsgs,lwork) 522 if ierror != 0: 523 msg = 'In return from call to shsgs in Spharmt.spectogrd ierror = %d' % ierror 524 raise ValueError, msg 525 else: 526 lwork = nlat*(nlon*nt+max(3*n2,nlon)) 527 datagrid, ierror = _spherepack.shsgc(nlon,a,b,self.wshsgc,lwork) 528 if ierror != 0: 529 msg = 'In return from call to shsgc in Spharmt.spectogrd ierror = %d' % ierror 530 raise ValueError, msg 531 532 if nt==1: 533 return numpy.squeeze(datagrid) 534 else: 535 return datagrid
536
537 - def getvrtdivspec(self, ugrid, vgrid, ntrunc=None):
538 539 """ 540 compute spectral coefficients of vorticity and divergence given vector wind. 541 542 @param ugrid: rank 2 or 3 numpy float32 array containing grid of zonal 543 winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number 544 of grids to be transformed. If ugrid is rank 2, nt is assumed to be 1. 545 546 @param vgrid: rank 2 or 3 numpy float32 array containing grid of meridional 547 winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number 548 of grids to be transformed. Both ugrid and vgrid must have the same shape. 549 550 @keyword ntrunc: optional spectral truncation limit. 551 (default self.nlat-1) 552 553 @return: C{B{vrtspec, divspec}} - rank 1 or 2 numpy complex arrays 554 of vorticity and divergence spherical harmonic coefficients with shape 555 shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt). 556 """ 557 558 # make sure ugrid,vgrid are rank 2 or 3 and same shape. 559 560 shapeu = ugrid.shape 561 shapev = vgrid.shape 562 563 if ntrunc is None: 564 ntrunc = self.nlat-1 565 566 if shapeu != shapev: 567 msg = 'getvrtdivspec input arrays must be same shape!' 568 raise ValueError, msg 569 570 571 if len(shapeu) !=2 and len(shapeu) !=3: 572 msg = 'getvrtdivspec needs rank two or three arrays!' 573 raise ValueError, msg 574 575 if shapeu[0] != self.nlat or shapeu[1] != self.nlon: 576 msg = 'getvrtdivspec needs input arrays whose first two dimensions are si%d and %d, got %d and %d' % (self.nlat, self.nlon, ugrid.shape[0], ugrid.shape[1],) 577 raise ValueError, msg 578 579 580 # check ntrunc. 581 582 if ntrunc < 0 or ntrunc+1 > shapeu[0]: 583 msg = 'ntrunc must be between 0 and %d' % (ugrid.shape[0]-1,) 584 raise ValueError, msg 585 586 nlat = self.nlat 587 nlon = self.nlon 588 if nlat%2: # nlat is odd 589 n2 = (nlat + 1)/2 590 else: 591 n2 = nlat/2 592 rsphere= self.rsphere 593 594 # convert from geographical to math coordinates, add extra dimension 595 # if necessary. 596 597 if len(shapeu) == 2: 598 nt = 1 599 w = numpy.reshape(ugrid, (nlat,nlon,1)) 600 v = -numpy.reshape(vgrid, (nlat,nlon,1)) 601 else: 602 nt = shapeu[2] 603 w = ugrid 604 v = -vgrid 605 606 # regular grid. 607 608 if self.gridtype == 'regular': 609 610 # vector harmonic analysis. 611 612 if self.legfunc == 'stored': 613 lwork = (2*nt+1)*nlat*nlon 614 br,bi,cr,ci,ierror = _spherepack.vhaes(v,w,self.wvhaes,lwork) 615 if ierror != 0: 616 msg = 'In return from call to vhaes in Spharmt.getvrtdivspec ierror = %d' % ierror 617 raise ValueError, msg 618 else: 619 lwork = nlat*(2*nt*nlon+max(6*n2,nlon)) 620 br,bi,cr,ci,ierror = _spherepack.vhaec(v,w,self.wvhaec,lwork) 621 if ierror != 0: 622 msg = 'In return from call to vhaec in Spharmt.getvrtdivspec ierror = %d' % ierror 623 raise ValueError, msg 624 625 # gaussian grid. 626 627 elif self.gridtype == 'gaussian': 628 629 # vector harmonic analysis. 630 631 if self.legfunc == 'stored': 632 lwork = (2*nt+1)*nlat*nlon 633 br,bi,cr,ci,ierror = _spherepack.vhags(v,w,self.wvhags,lwork) 634 if ierror != 0: 635 msg = 'In return from call to vhags in Spharmt.getvrtdivspec ierror = %d' % ierror 636 raise ValueError, msg 637 else: 638 lwork = 2*nlat*(2*nlon*nt+3*n2) 639 br,bi,cr,ci,ierror = _spherepack.vhagc(v,w,self.wvhagc,lwork) 640 if ierror != 0: 641 msg = 'In return from call to vhagc in Spharmt.getvrtdivspec ierror = %d' % ierror 642 raise ValueError, msg 643 644 # convert vector harmonic coeffs to 1d complex coefficients 645 # of vorticity and divergence. 646 647 vrtspec, divspec = _spherepack.twodtooned_vrtdiv(br,bi,cr,ci,ntrunc,rsphere) 648 649 if nt==1: 650 return numpy.squeeze(vrtspec), numpy.squeeze(divspec) 651 else: 652 return vrtspec, divspec
653
654 - def getuv(self, vrtspec, divspec):
655 656 """ 657 compute vector wind on grid given complex spectral coefficients 658 of vorticity and divergence. 659 660 @param vrtspec: rank 1 or 2 numpy complex array of vorticity spectral 661 coefficients, with shape (ntrunc+1)*(ntrunc+2)/2 or 662 ((ntrunc+1)*(ntrunc+2)/2,nt) (where ntrunc is the triangular truncation 663 and nt is the number of spectral arrays to be transformed). 664 If vrtspec is rank 1, nt is assumed to be 1. 665 666 @param divspec: rank 1 or 2 numpy complex array of divergence spectral 667 coefficients, with shape (ntrunc+1)*(ntrunc+2)/2 or 668 ((ntrunc+1)*(ntrunc+2)/2,nt) (where ntrunc is the triangular truncation 669 and nt is the number of spectral arrays to be transformed). 670 Both vrtspec and divspec must have the same shape. 671 672 @return: C{B{ugrid, vgrid}} - rank 2 or 3 numpy float32 arrays containing 673 gridded zonal and meridional winds. Shapes are either (nlat,nlon) or 674 (nlat,nlon,nt). 675 """ 676 677 shapevrt = vrtspec.shape 678 shapediv = divspec.shape 679 680 # make sure vrtspec, divspec are rank 1 or 2, and have the same shape. 681 682 if shapevrt != shapediv: 683 msg = 'vrtspec, divspec must be same size in getuv!' 684 raise ValueError, msg 685 686 if len(shapevrt) !=1 and len(shapevrt) !=2: 687 msg = 'getuv needs rank one or two input arrays!' 688 raise ValueError, msg 689 690 # infer ntrunc from size of dataspec (dataspec must be rank 1!) 691 # dataspec is assumed to have size (ntrunc+1)*(ntrunc+2)/2 692 693 ntrunc = int(-1.5 + 0.5*math.sqrt(9.-8.*(1.-vrtspec.shape[0]))) 694 695 nlat = self.nlat 696 nlon = self.nlon 697 if nlat%2: # nlat is odd 698 n2 = (nlat + 1)/2 699 else: 700 n2 = nlat/2 701 rsphere= self.rsphere 702 703 if len(vrtspec.shape) == 1: 704 nt = 1 705 vrtspec = numpy.reshape(vrtspec, (vrtspec.shape[0],1)) 706 divspec = numpy.reshape(divspec, (divspec.shape[0],1)) 707 else: 708 nt = vrtspec.shape[1] 709 710 # convert 1d complex arrays of vort, div to 2d vector harmonic arrays. 711 712 br,bi,cr,ci = _spherepack.onedtotwod_vrtdiv(vrtspec,divspec,nlat,rsphere) 713 714 # regular grid. 715 716 if self.gridtype == 'regular': 717 718 # vector harmonic synthesis. 719 720 if self.legfunc == 'stored': 721 lwork = (2*nt+1)*nlat*nlon 722 v, w, ierror = _spherepack.vhses(nlon,br,bi,cr,ci,self.wvhses,lwork) 723 if ierror != 0: 724 msg = 'In return from call to vhses in Spharmt.getuv ierror = %d' % ierror 725 raise ValueError, msg 726 else: 727 lwork = nlat*(2*nt*nlon+max(6*n2,nlon)) 728 v, w, ierror = _spherepack.vhsec(nlon,br,bi,cr,ci,self.wvhsec,lwork) 729 if ierror != 0: 730 msg = 'In return from call to vhsec in Spharmt.getuv ierror = %d' % ierror 731 raise ValueError, msg 732 733 # gaussian grid. 734 735 elif self.gridtype == 'gaussian': 736 737 # vector harmonic synthesis. 738 739 if self.legfunc == 'stored': 740 lwork = (2*nt+1)*nlat*nlon 741 v, w, ierror = _spherepack.vhsgs(nlon,br,bi,cr,ci,self.wvhsgs,lwork) 742 if ierror != 0: 743 msg = 'In return from call to vhsgs in Spharmt.getuv ierror = %d' % ierror 744 raise ValueError, msg 745 else: 746 lwork = nlat*(2*nt*nlon+max(6*n2,nlon)) 747 v, w, ierror = _spherepack.vhsgc(nlon,br,bi,cr,ci,self.wvhsgc,lwork) 748 if ierror != 0: 749 msg = 'In return from call to vhsgc in Spharmt.getuv ierror = %d' % ierror 750 raise ValueError, msg 751 752 # convert to u and v in geographical coordinates. 753 754 if nt == 1: 755 return numpy.reshape(w, (nlat,nlon)), -numpy.reshape(v, (nlat,nlon)) 756 else: 757 return w,-v
758
759 - def getpsichi(self, ugrid, vgrid, ntrunc=None):
760 761 """ 762 compute streamfunction and velocity potential on grid given vector wind. 763 764 @param ugrid: rank 2 or 3 numpy float32 array containing grid of zonal 765 winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number 766 of grids to be transformed. If ugrid is rank 2, nt is assumed to be 1. 767 768 @param vgrid: rank 2 or 3 numpy float32 array containing grid of meridional 769 winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number 770 of grids to be transformed. Both ugrid and vgrid must have the same shape. 771 772 @keyword ntrunc: optional spectral truncation limit. 773 (default self.nlat-1) 774 775 @return: C{B{psigrid, chigrid}} - rank 2 or 3 numpy float32 arrays 776 of gridded streamfunction and velocity potential. Shapes are either 777 (nlat,nlon) or (nlat,nlon,nt). 778 """ 779 780 # make sure ugrid,vgrid are rank 2 or 3 and same shape. 781 782 shapeu = ugrid.shape 783 shapev = vgrid.shape 784 785 if ntrunc is None: 786 ntrunc = self.nlat-1 787 788 if shapeu != shapev: 789 msg = 'getvrtdivspec input arrays must be same shape!' 790 raise ValueError, msg 791 792 793 if len(shapeu) !=2 and len(shapeu) !=3: 794 msg = 'getvrtdivspec needs rank two or three arrays!' 795 raise ValueError, msg 796 797 if shapeu[0] != self.nlat or shapeu[1] != self.nlon: 798 msg = 'getpsichi needs input arrays whose first two dimensions are si%d and %d, got %d and %d' % (self.nlat, self.nlon, ugrid.shape[0], ugrid.shape[1],) 799 raise ValueError, msg 800 801 # check ntrunc. 802 803 if ntrunc < 0 or ntrunc+1 > ugrid.shape[0]: 804 msg = 'ntrunc must be between 0 and %d' % (ugrid.shape[0]-1,) 805 raise ValueError, msg 806 807 # compute spectral coeffs of vort, div. 808 809 vrtspec, divspec = self.getvrtdivspec(ugrid, vgrid, ntrunc) 810 811 # number of grids to compute. 812 813 if len(vrtspec.shape) == 1: 814 nt = 1 815 vrtspec = numpy.reshape(vrtspec, ((ntrunc+1)*(ntrunc+2)/2,1)) 816 divspec = numpy.reshape(divspec, ((ntrunc+1)*(ntrunc+2)/2,1)) 817 else: 818 nt = vrtspec.shape[1] 819 820 # convert to spectral coeffs of psi, chi. 821 822 psispec = _spherepack.invlap(vrtspec, self.rsphere) 823 chispec = _spherepack.invlap(divspec, self.rsphere) 824 825 # inverse transform to grid. 826 827 psigrid = self.spectogrd(psispec) 828 chigrid = self.spectogrd(chispec) 829 830 return psigrid, chigrid
831
832 - def getgrad(self, chispec):
833 834 """ 835 compute vector gradient on grid given complex spectral coefficients. 836 837 @param chispec: rank 1 or 2 numpy complex array with shape 838 (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing 839 complex spherical harmonic coefficients (where ntrunc is the 840 triangular truncation limit and nt is the number of spectral arrays 841 to be transformed). If chispec is rank 1, nt is assumed to be 1. 842 843 @return: C{B{uchi, vchi}} - rank 2 or 3 numpy float32 arrays containing 844 gridded zonal and meridional components of the vector gradient. 845 Shapes are either (nlat,nlon) or (nlat,nlon,nt). 846 """ 847 848 # make sure chispec is rank 1 or 2. 849 850 if len(chispec.shape) !=1 and len(chispec.shape) !=2: 851 msg = 'getgrad needs rank one or two arrays!' 852 raise ValueError, msg 853 854 # infer ntrunc from size of chispec (chispec must be rank 1!) 855 # chispec is assumed to have size (ntrunc+1)*(ntrunc+2)/2 856 857 ntrunc = int(-1.5 + 0.5*math.sqrt(9.-8.*(1.-chispec.shape[0]))) 858 859 # number of grids to compute. 860 861 if len(chispec.shape) == 1: 862 nt = 1 863 chispec = numpy.reshape(chispec, ((ntrunc+1)*(ntrunc+2)/2,1)) 864 else: 865 nt = chispec.shape[1] 866 867 # convert chispec to divspec. 868 869 divspec = _spherepack.lap(chispec,self.rsphere) 870 871 # call getuv, with vrtspec=0, to get uchi,vchi. 872 873 chispec[:,:] = 0. 874 875 uchi, vchi = self.getuv(chispec, divspec) 876 877 return uchi, vchi
878
879 - def specsmooth(self, datagrid, smooth):
880 881 """ 882 isotropic spectral smoothing on a sphere. 883 884 @param datagrid: rank 2 or 3 numpy float32 array with shape (nlat,nlon) or 885 (nlat,nlon,nt), where nt is the number of grids to be smoothed. If 886 datagrid is rank 2, nt is assumed to be 1. 887 888 @param smooth: rank 1 array of length nlat containing smoothing factors 889 as a function of total wavenumber. 890 891 @return: C{B{datagrid}} - rank 2 or 3 numpy float32 array with shape 892 (nlat,nlon) or (nlat,nlon,nt) containing the smoothed grids. 893 """ 894 895 # check that datagrid is rank 2 or 3 with size (self.nlat, self.nlon) or 896 # (self.nlat, self.nlon, nt) where nt is number of grids to transform. 897 898 if len(datagrid.shape) > 3: 899 msg = 'specsmooth needs a rank two or three array, got %d' % (len(datagrid.shape),) 900 raise ValueError, msg 901 902 if datagrid.shape[0] != self.nlat or datagrid.shape[1] != self.nlon: 903 msg = 'specsmooth needs an array of size %d by %d, got %d by %d' % (self.nlat, self.nlon, datagrid.shape[0], datagrid.shape[1],) 904 raise ValueError, msg 905 906 907 # make sure smooth is rank 1, same size as datagrid.shape[0] 908 909 if len(smooth.shape) !=1 or smooth.shape[0] != datagrid.shape[0]: 910 msg = 'smooth must be rank 1 and same size as datagrid.shape[0] in specsmooth!' 911 raise ValueError, msg 912 913 # grid to spectral transform. 914 915 nlat = self.nlat 916 dataspec = self.grdtospec(datagrid, nlat-1) 917 918 # multiply spectral coeffs. by smoothing factor. 919 920 dataspec = _spherepack.multsmoothfact(dataspec, smooth) 921 922 # spectral to grid transform. 923 924 datagrid = self.spectogrd(dataspec) 925 926 return datagrid
927
928 -def regrid(grdin, grdout, datagrid, ntrunc=None, smooth=None):
929 """ 930 regrid data using spectral interpolation, while performing 931 optional spectral smoothing and/or truncation. 932 933 @param grdin: Spharmt class instance describing input grid. 934 935 @param grdout: Spharmt class instance describing output grid. 936 937 @param datagrid: data on input grid (grdin.nlat x grdin.nlon). If 938 datagrid is rank 3, last dimension is the number of grids to interpolate. 939 940 @keyword ntrunc: optional spectral truncation limit for datagrid 941 (default min(grdin.nlat-1,grdout.nlat-1)). 942 943 @keyword smooth: rank 1 array of length grdout.nlat containing smoothing 944 factors as a function of total wavenumber (default is no smoothing). 945 946 @return: C{B{datagrid}} - interpolated (and optionally smoothed) array(s) 947 on grdout.nlon x grdout.nlat grid. 948 """ 949 950 # check that datagrid is rank 2 or 3 with size (grdin.nlat, grdin.nlon) or 951 # (grdin.nlat, grdin.nlon, nt) where nt is number of grids to transform. 952 953 if len(datagrid.shape) > 3: 954 msg = 'regrid needs a rank two or three array, got %d' % (len(datagrid.shape),) 955 raise ValueError, msg 956 957 if datagrid.shape[0] != grdin.nlat or datagrid.shape[1] != grdin.nlon: 958 msg = 'grdtospec needs an array of size %d by %d, got %d by %d' % (grdin.nlat, grdin.nlon, datagrid.shape[0], datagrid.shape[1],) 959 raise ValueError, msg 960 961 if smooth is not None and (len(smooth.shape) !=1 or smooth.shape[0] != grdout.nlat): 962 msg = 'smooth must be rank 1 size grdout.nlat in regrid!' 963 raise ValueError, msg 964 965 if ntrunc is None: 966 ntrunc = min(grdout.nlat-1,grdin.nlat-1) 967 968 dataspec = grdin.grdtospec(datagrid,ntrunc) 969 970 if smooth is not None: 971 dataspec = _spherepack.multsmoothfact(dataspec, smooth) 972 973 return grdout.spectogrd(dataspec)
974
975 -def gaussian_lats_wts(nlat):
976 977 """ 978 compute the gaussian latitudes (in degrees) and quadrature weights. 979 980 @param nlat: number of gaussian latitudes desired. 981 982 @return: C{B{lats, wts}} - rank 1 numpy float64 arrays containing 983 gaussian latitudes (in degrees north) and gaussian quadrature weights. 984 """ 985 986 # get the gaussian colatitudes and weights using gaqd. 987 988 colats, wts, ierror = _spherepack.gaqd(nlat) 989 990 if ierror != 0: 991 msg = 'In return from call to gaqd ierror = %d' % ierror 992 raise ValueError, msg 993 994 # convert to degrees north latitude. 995 996 lats = 90.0 - colats*180.0/math.pi 997 998 return lats, wts
999
1000 -def getspecindx(ntrunc):
1001 1002 """ 1003 compute indices of zonal wavenumber (indxm) and degree (indxn) 1004 for complex spherical harmonic coefficients. 1005 1006 @param ntrunc: spherical harmonic triangular truncation limit. 1007 1008 @return: C{B{indxm, indxn}} - rank 1 numpy Int32 arrays 1009 containing zonal wavenumber (indxm) and degree (indxn) of 1010 spherical harmonic coefficients. 1011 """ 1012 1013 indexn = numpy.indices((ntrunc+1,ntrunc+1))[1,:,:] 1014 indexm = numpy.indices((ntrunc+1,ntrunc+1))[0,:,:] 1015 indices = numpy.nonzero(numpy.greater(indexn, indexm-1).flatten()) 1016 indxn = numpy.take(indexn.flatten(),indices) 1017 indxm = numpy.take(indexm.flatten(),indices) 1018 1019 return numpy.squeeze(indxm), numpy.squeeze(indxn)
1020
1021 -def getgeodesicpts(m):
1022 """ 1023 computes the lat/lon values of the points on the surface of the sphere 1024 corresponding to a twenty-sided (icosahedral) geodesic. 1025 1026 @param m: the number of points on the edge of a single geodesic triangle. 1027 There are 10*(m-1)**2+2 total geodesic points, including the poles. 1028 1029 @return: C{B{lats, lons}} - rank 1 numpy float32 arrays containing 1030 the latitudes and longitudes of the geodesic points (in degrees). These 1031 points are nearly evenly distributed on the surface of the sphere. 1032 """ 1033 x,y,z = _spherepack.ihgeod(m) 1034 # convert cartesian coords to lat/lon. 1035 rad2dg = 180./math.pi 1036 r1 = x*x+y*y 1037 r = numpy.sqrt(r1+z*z) 1038 r1 = numpy.sqrt(r1) 1039 xtmp = numpy.where(numpy.logical_or(x,y),x,numpy.ones(x.shape,numpy.float32)) 1040 ztmp = numpy.where(numpy.logical_or(r1,z),z,numpy.ones(z.shape,numpy.float32)) 1041 lons = rad2dg*numpy.arctan2(y,xtmp)+180. 1042 lats = rad2dg*numpy.arctan2(r1,ztmp)-90. 1043 lat = numpy.zeros(10*(m-1)**2+2,numpy.float32) 1044 lon = numpy.zeros(10*(m-1)**2+2,numpy.float32) 1045 # first two points are poles. 1046 lat[0] = 90; lat[1] = -90. 1047 lon[0] = 0.; lon[1] = 0. 1048 lat[2:] = lats[0:2*(m-1),0:m-1,:].flatten() 1049 lon[2:] = lons[0:2*(m-1),0:m-1,:].flatten() 1050 return lat,lon
1051
1052 -def legendre(lat,ntrunc):
1053 """ 1054 calculate associated legendre functions for triangular truncation T(ntrunc), 1055 at a given latitude. 1056 1057 @param lat: the latitude (in degrees) to compute the associate legendre 1058 functions. 1059 1060 @param ntrunc: the triangular truncation limit. 1061 1062 @return: C{B{pnm}} - rank 1 numpy float32 array containing the 1063 (C{B{ntrunc}}+1)*(C{B{ntrunc}}+2)/2 associated legendre functions at 1064 latitude C{B{lat}}. 1065 """ 1066 return _spherepack.getlegfunc(lat,ntrunc)
1067
1068 -def specintrp(lon,dataspec,legfuncs):
1069 """ 1070 spectral interpolation given spherical harmonic coefficients. 1071 1072 @param lon: longitude (in degrees) of point on a sphere to interpolate to. 1073 1074 @param dataspec: spectral coefficients of function to interpolate. 1075 1076 @param legfuncs: associated legendre functions with same triangular 1077 truncation as C{B{dataspec}} (computed using L{legendre}), computed 1078 at latitude of interpolation point. 1079 1080 @return: C{B{ob}} - interpolated value. 1081 """ 1082 ntrunc1 = int(-1.5 + 0.5*math.sqrt(9.-8.*(1.-dataspec.shape[0]))) 1083 ntrunc2 = int(-1.5 + 0.5*math.sqrt(9.-8.*(1.-legfuncs.shape[0]))) 1084 if ntrunc1 != ntrunc2: 1085 raise ValueError, 'first dimensions of dataspec and legfuncs in Spharmt.specintrp imply inconsistent spectral truncations - they must be the same!' 1086 return _spherepack.specintrp((math.pi/180.)*lon,ntrunc1,dataspec,legfuncs)
1087