MPAS-JEDI
binning_utils.py
Go to the documentation of this file.
1 #!/usr/bin/env python3
2 
3 from collections import defaultdict
4 from collections.abc import Iterable
5 from copy import deepcopy
6 import inspect
7 import logging
8 import numpy as np
9 import os
10 import plot_utils as pu
11 from binning_params import allSCIErrParams, ABEIParams
12 import var_utils as vu
13 
14 _logger = logging.getLogger(__name__)
15 
16 #===========================
17 # binning method descriptors
18 #===========================
19 # TODO(JJG): should be moved to binning_params
20 
21 # identity
22 identityBinMethod = 'identity'
23 
24 # none
25 noBinMethod = identityBinMethod
26 
27 #QC
28 goodQCMethod = 'good'
29 badQCMethod = 'bad'
30 allQCMethod = 'all'
31 
32 #jet-stream pressure
33 P_jet_min = 250.0
34 P_jet_max = 350.0
35 P_jet_val = '{:.0f}'.format(0.5 * (P_jet_min + P_jet_max))
36 PjetMethod = 'P='+P_jet_val+'hPa'
37 
38 #jet-stream altitude
39 alt_jet_min = 9500.0
40 alt_jet_max = 10500.0
41 alt_jet_val = '{:.0f}'.format(0.5 * (alt_jet_min + alt_jet_max))
42 altjetMethod = 'alt='+alt_jet_val+'m'
43 
44 #LocalHour
45 LH0 = 0.0
46 LH1 = 23.0
47 LHDT = 1.0
48 
49 #named latitude-bands
50 latbandsMethod = 'LatBands'
51 
52 #named surface-type-bands
53 seasurfMethod = 'sea'
54 landsurfMethod = 'land'
55 mixsurfMethod = 'mixed-land-sea'
56 allsurfMethod = 'all-surface'
57 surfbandsMethod = 'surface-type'
58 
59 #named cloudiness-bands
60 clrskyMethod = 'clear'
61 cldskyMethod = 'cloudy'
62 mixskyMethod = 'mixed-clrcld'
63 allskyMethod = 'allsky'
64 cloudbandsMethod = 'cloudiness'
65 
66 # symmetric cloud impact (SCI)
67 OkamotoMethod = 'Okamoto'
68 ScaleOkamotoMethod = 'ScaledOkamoto'
69 ModHarnischMethod = 'ModHarnisch'
70 ScaleModHarnischMethod = 'ScaledModHarnisch'
71 
72 # geostationary IR instrument longitude parameters
73 geoirlatlonboxMethod = 'geoirBox'
74 geoirMaxZenith = 65.0
75 
76 geoirlatlonBoxParams = defaultdict(list)
77 
78 # activate desired instruments by uncommenting them below
79 
80 #seviri_m11 = 'seviri_m11'
81 #geoirlatlonBoxParams['values'] += [seviri_m11]
82 #geoirlatlonBoxParams['centerLon'] += [0.]
83 
84 #seviri_m08 = 'seviri_m08'
85 #geoirlatlonBoxParams['values'] += [seviri_m08]
86 #geoirlatlonBoxParams['centerLon'] += [41.5]
87 
88 ahi_himawari8 = 'ahi_himawari8'
89 geoirlatlonBoxParams['values'] += [ahi_himawari8]
90 geoirlatlonBoxParams['centerLon'] += [140.7]
91 
92 #abi_g17 = 'abi_g17'
93 #geoirlatlonBoxParams['values'] += [abi_g17]
94 #geoirlatlonBoxParams['centerLon'] += [360. - 137.2]
95 
96 abi_g16 = 'abi_g16'
97 geoirlatlonBoxParams['values'] += [abi_g16]
98 geoirlatlonBoxParams['centerLon'] += [360. - 75.2]
99 
100 # glint angle
101 maxGlint = 90.0
102 
103 #========================
104 # binning where functions
105 #========================
106 
107 # Note: NaN and Inf values have mask set to true by default
108 
109 def equalBound(x, bound, missingValue=True):
110  mask = np.empty_like(x, bool)
111  if isinstance(bound, str):
112  mask = np.char.equal(x, bound)
113  else:
114  finite = np.isfinite(x)
115  mask[finite] = np.equal(x[finite], bound)
116  mask[~finite] = missingValue
117  return mask
118 
119 def notEqualBound(x, bound, missingValue=True):
120  mask = np.empty_like(x, bool)
121  if isinstance(bound, str):
122  mask = np.char.not_equal(x, bound)
123  else:
124  finite = np.isfinite(x)
125  mask[finite] = np.not_equal(x[finite], bound)
126  mask[~finite] = missingValue
127  return mask
128 
129 def notEqualAnyBound(x, bounds, missingValue=True):
130  assert isinstance(bounds, Iterable), \
131  ('ERROR, bounds must be Iterable for notEqualAnyBound')
132 
133  mask = np.full_like(x, True, bool)
134  if isinstance(bounds[0], str):
135  for bound in bounds:
136  mask = np.logical_and(mask,
137  np.char.not_equal(x, bound))
138  else:
139  finite = np.isfinite(x)
140  for bound in bounds:
141  mask[finite] = np.logical_and(mask[finite],
142  np.not_equal(x[finite], bound))
143  mask[~finite] = missingValue
144  return mask
145 
146 def lessEqualBound(x, bound, missingValue=True):
147  finite = np.isfinite(x)
148  mask = np.empty_like(finite, bool)
149  mask[finite] = np.less_equal(x[finite], bound)
150  mask[~finite] = missingValue
151  return mask
152 
153 def lessBound(x, bound, missingValue=True):
154  finite = np.isfinite(x)
155  mask = np.empty_like(finite, bool)
156  mask[finite] = np.less(x[finite], bound)
157  mask[~finite] = missingValue
158  return mask
159 
160 def greatEqualBound(x, bound, missingValue=True):
161  finite = np.isfinite(x)
162  mask = np.empty_like(finite, bool)
163  mask[finite] = np.greater_equal(x[finite], bound)
164  mask[~finite] = missingValue
165  return mask
166 
167 def greatBound(x, bound, missingValue=True):
168  finite = np.isfinite(x)
169  mask = np.empty_like(finite, bool)
170  mask[finite] = np.greater(x[finite], bound)
171  mask[~finite] = missingValue
172  return mask
173 
174 def betweenBounds(x, bound1, bound2, missingValue=True):
175  finite = np.isfinite(x)
176  belowbounds = lessEqualBound(x, bound1)
177  abovebounds = greatEqualBound(x, bound2)
178  mask = np.logical_not(np.logical_or(
179  belowbounds, abovebounds ))
180  mask[~finite] = missingValue
181  return mask
182 
183 
184 #=========================================================
185 # ObsFunction classes to be accessed w/ ObsFunctionWrapper
186 # i.e., functions of variables contained in the database
187 #=========================================================
188 
190  def __init__(self):
191  self.baseVarsbaseVars = []
192  self.baseVarsbaseVars.append(vu.senzenMeta)
193  self.baseVarsbaseVars.append(vu.senaziMeta)
194  self.baseVarsbaseVars.append(vu.solzenMeta)
195  self.baseVarsbaseVars.append(vu.solaziMeta)
196 
197  def evaluate(self, dbVals, insituParameters):
198  senazi = dbVals[insituParameters[vu.senaziMeta]]
199  solazi = dbVals[insituParameters[vu.solaziMeta]]
200 
201  relazi = np.abs(np.subtract(solazi,senazi))
202  p = greatBound(relazi, 180.0, False)
203  relazi[p] = np.subtract(360.0, relazi[p])
204  relazi = np.multiply(np.subtract(180.0, relazi), vu.deg2rad)
205 
206  senzen = np.multiply(dbVals[insituParameters[vu.senzenMeta]], vu.deg2rad)
207  solzen = np.multiply(dbVals[insituParameters[vu.solzenMeta]], vu.deg2rad)
208 
209  glint = np.add(np.multiply(np.cos(solzen), np.cos(senzen)),
210  np.multiply(np.sin(solzen),
211  np.multiply(np.sin(senzen), np.cos(relazi))))
212 
213  glint[greatBound(glint, 1.0)] = np.NaN
214  glint[lessBound(glint, -1.0)] = np.NaN
215 
216  glint = np.multiply(np.arccos(glint), vu.rad2deg)
217  glint[greatBound(glint, maxGlint, False)] = maxGlint
218 
219  return glint
220 
221 
222 class LocalHour:
223  def __init__(self):
224  self.baseVarsbaseVars = []
225  self.baseVarsbaseVars.append(vu.datetimeMeta)
226  self.baseVarsbaseVars.append(vu.lonMeta)
227 
228  def evaluate(self, dbVals, insituParameters):
229  TimeStr = dbVals[insituParameters[vu.datetimeMeta]]
230  tzOffset = np.divide(dbVals[insituParameters[vu.lonMeta]], 15.0)
231 
232  hh = np.empty_like(tzOffset, dtype=np.float32)
233  mmi = np.empty_like(tzOffset, dtype=np.float32)
234  ss = np.empty_like(tzOffset, dtype=np.float32)
235 
236  t0 = LH0
237  t1 = LH1
238  dt = LHDT
239  for ii, Time in enumerate(TimeStr):
240  ## Expecting Time to fit YYYY-MM-DDThh:mm:ssZ
241  # YYYY[ii] = float(Time[0:4])
242  # MMo[ii] = float(Time[5:7])
243  # DD[ii] = float(Time[8:10])
244  hh[ii] = float(Time[11:13])
245  mmi[ii] = float(Time[14:16]) / 60.0
246  ss[ii] = float(Time[17:19]) / 3600.0
247 
248  LH = hh + mmi + ss + tzOffset
249 
250  yesterday = (LH < t0 - 0.5*dt)
251  LH[yesterday] = LH[yesterday] + 24.0
252 
253  LH = np.mod(LH, 24.0)
254 
255  tomorrow = (LH >= t1 + 0.5*dt)
256  LH[tomorrow] = LH[tomorrow] - 24.0
257 
258  return LH
259 
260 
261 #classes related to the Asymmetric Cloud Impact (ACI)
263  def __init__(self):
264  self.baseVarsbaseVars = []
265  self.baseVarsbaseVars.append(vu.selfObsValue)
266  # self.baseVars.append(vu.selfDepValue)
267  self.baseVarsbaseVars.append(vu.selfHofXValue)
268  self.baseVarsbaseVars.append(vu.clrskyBTDiag)
269 
270  def evaluate(self,dbVals,insituParameters):
271  # Minamide and Zhang (2018)
272  BTobs = dbVals[insituParameters[vu.selfObsValue]]
273  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
274  # BTbak = np.add(BTdep,BTobs)
275  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
276  BTclr = deepcopy(dbVals[insituParameters[vu.clrskyBTDiag]])
277  p = lessBound(BTclr, 1.0, False)
278  BTclr[p] = BTbak[p]
279  # Minamide and Zhange (2018), maybe run into problems with non-clear-sky-bias-corrected
280  ACI = np.subtract(np.abs(np.subtract(BTobs, BTclr)),
281  np.abs(np.subtract(BTbak, BTclr)))
282 # some ideas, but still causes bias in ACI for non-clear-sky-bias-corrected
283 # zeros = np.zeros(BTbak.shape)
284 # ACI = np.subtract(np.maximum(zeros, np.subtract(BTclr, BTobs)),
285 # np.abs(np.subtract(BTbak, BTclr)))
286 ## np.maximum(zeros, np.subtract(BTclr, BTbak)))
287 
288  return ACI
289 
290 
292  minLambda = 1.0
293  maxLambda = 1.4
294  def __init__(self):
296  self.baseVarsbaseVars = []
297  self.baseVarsbaseVars = pu.uniqueMembers(self.baseVarsbaseVars + self.ACIACI.baseVars)
298 
299  def evaluate(self, dbVals, insituParameters):
300  osName = insituParameters['osName']
301  if osName is None or osName not in ABEIParams:
302  _logger.error('osName not available in ABEIParams => '+osName)
303  os._exit(1)
304 
305  # varName, ch = vu.splitIntSuffix(insituParameters[vu.selfDepValue])
306  varName, ch = vu.splitIntSuffix(insituParameters[vu.selfHofXValue])
307  LambdaOverACI = ABEIParams[osName][(int(ch))]['LambdaOverACI']
308 
309  ACI = self.ACIACI.evaluate(dbVals, insituParameters)
310  lambdaOut = np.ones(ACI.shape)
311  crit = (ACI > 0.0)
312  lambdaOut[crit] = np.multiply(ACI[crit], LambdaOverACI) + self.minLambdaminLambda
313  crit = (ACI >= (self.maxLambdamaxLambda - self.minLambdaminLambda) / LambdaOverACI)
314  lambdaOut[crit] = self.maxLambdamaxLambda
315 
316  return lambdaOut
317 
318 
319 #classes related to the Symmetric Cloud Impact (SCI)
321  def __init__(self):
322  self.baseVarsbaseVars = []
323  self.baseVarsbaseVars.append(vu.selfObsValue)
324  # self.baseVars.append(vu.selfDepValue)
325  self.baseVarsbaseVars.append(vu.selfHofXValue)
326  self.baseVarsbaseVars.append(vu.clrskyBTDiag)
327 
328  def evaluate(self, dbVals, insituParameters):
329  # Okamoto, et al.
330  # Co = abs(Bias-Corrected BTobs - BTclr)
331  # Cm = abs(BTbak - BTclr)
332 
333  BTobs = dbVals[insituParameters[vu.selfObsValue]]
334  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
335  # BTbak = np.add(BTdep, BTobs)
336  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
337  BTclr = deepcopy(dbVals[insituParameters[vu.clrskyBTDiag]])
338  BTclr[BTclr < 1.0] = BTbak[BTclr < 1.0]
339  SCI = np.multiply( 0.5,
340  np.add(np.abs(np.subtract(BTobs, BTclr)),
341  np.abs(np.subtract(BTbak, BTclr))) )
342  return SCI
343 
344 
346  def __init__(self):
347  self.baseVarsbaseVars = []
348  self.baseVarsbaseVars.append(vu.selfObsValue)
349  # self.baseVars.append(vu.selfDepValue)
350  self.baseVarsbaseVars.append(vu.selfHofXValue)
351  self.baseVarsbaseVars.append(vu.clrskyBTDiag)
352  self.baseVarsbaseVars.append(vu.cldfracMeta)
353 
354  def evaluate(self, dbVals, insituParameters):
355  BTobs = dbVals[insituParameters[vu.selfObsValue]]
356  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
357  # BTbak = np.add(BTdep, BTobs)
358  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
359  BTclr = deepcopy(dbVals[insituParameters[vu.clrskyBTDiag]])
360  BTclr[BTclr < 1.0] = BTbak[BTclr < 1.0]
361  CldFrac = dbVals[insituParameters[vu.cldfracMeta]]
362 
363  #Scale both Co and Cm by retrieved cloud fraction
364  # SCI = np.multiply( 0.5,
365  # np.multiply(CldFrac,
366  # np.add(np.abs(np.subtract(BTobs, BTclr)),
367  # np.abs(np.subtract(BTbak, BTclr))) ) )
368 
369  #Scale only Co by retrieved cloud fraction
370  SCI = np.multiply( 0.5,
371  np.add(
372  np.multiply(CldFrac,
373  np.abs(np.subtract(BTobs, BTclr))),
374  np.abs(np.subtract(BTbak, BTclr)) ) )
375 
376  return SCI
377 
378 
380  def __init__(self):
381  self.baseVarsbaseVars = []
382  self.baseVarsbaseVars.append(vu.selfObsValue)
383  # self.baseVars.append(vu.selfDepValue)
384  self.baseVarsbaseVars.append(vu.selfHofXValue)
385  self.baseVarsbaseVars.append(vu.clrskyBTDiag)
386 
387  def evaluate(self, dbVals, insituParameters):
388  # Modified Harnisch, et al.
389  BTobs = dbVals[insituParameters[vu.selfObsValue]]
390  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
391  # BTbak = np.add(BTdep, BTobs)
392  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
393  BTclr = deepcopy(dbVals[insituParameters[vu.clrskyBTDiag]])
394  BTclr[BTclr < 1.0] = BTbak[BTclr < 1.0]
395  zeros = np.full_like(BTbak,0.0)
396  SCI = np.multiply( 0.5,
397  np.add(np.maximum(zeros, np.subtract(BTclr, BTobs)),
398  np.maximum(zeros, np.subtract(BTclr, BTbak))) )
399  return SCI
400 
401 
403  def __init__(self):
404  self.baseVarsbaseVars = []
405  self.baseVarsbaseVars.append(vu.selfObsValue)
406  # self.baseVars.append(vu.selfDepValue)
407  self.baseVarsbaseVars.append(vu.selfHofXValue)
408  self.baseVarsbaseVars.append(vu.clrskyBTDiag)
409  self.baseVarsbaseVars.append(vu.cldfracMeta)
410 
411  def evaluate(self, dbVals, insituParameters):
412  # Modified Harnisch, et al.
413  BTobs = dbVals[insituParameters[vu.selfObsValue]]
414  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
415  # BTbak = np.add(BTdep, BTobs)
416  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
417  BTclr = deepcopy(dbVals[insituParameters[vu.clrskyBTDiag]])
418  BTclr[BTclr < 1.0] = BTbak[BTclr < 1.0]
419  CldFrac = dbVals[insituParameters[vu.cldfracMeta]]
420 
421  zeros = np.full_like(BTbak,0.0)
422  SCI = np.multiply( 0.5,
423  np.multiply(CldFrac,
424  np.add(np.maximum(zeros, np.subtract(BTclr, BTobs)),
425  np.maximum(zeros, np.subtract(BTclr, BTbak))) ) )
426  return SCI
427 
428 
430  def __init__(self):
431  self.baseVarsbaseVars = []
432  self.baseVarsbaseVars.append(vu.selfObsValue)
433  # self.baseVars.append(vu.selfDepValue)
434  self.baseVarsbaseVars.append(vu.selfHofXValue)
435  self.baseVarsbaseVars.append(vu.selfErrorValue)
436 
437  def evaluate(self, dbVals, insituParameters):
438  BTerr = dbVals[insituParameters[vu.selfErrorValue]]
439  BTerr[BTerr==0.0] = np.NaN
440  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
441  BTobs = dbVals[insituParameters[vu.selfObsValue]]
442  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
443  BTdep = np.subtract(BTbak, BTobs)
444 
445  return np.divide(BTdep, BTerr)
446 
447 
448 mpasFCRes = 120 # km [30, 120]
449 biasCorrectType = {} #[None, 'constant', 'varbc']
450 biasCorrectType['abi_g16'] = 'constant'
451 biasCorrectType['ahi_himawari8'] = None
452 
454  def __init__(self):
455  self.baseVarsbaseVars = []
456  # self.baseVars.append(vu.selfDepValue)
457  self.baseVarsbaseVars.append(vu.selfObsValue)
458  self.baseVarsbaseVars.append(vu.selfHofXValue)
459 
460  def evaluate(self, dbVals, insituParameters, SCISTDName, SCI):
461  # Parameterize BTerr as a ramped step function
462  # --------------------------------------------
463  # STD1 |. . . ____
464  # | /.
465  # | / .
466  # | / .
467  # STD0 |__/ .
468  # | . .
469  # |__.___.___
470  # . .
471  # SCI0 SCI1
472  #---------------------------------------------
473  osName = insituParameters['osName']
474  SCIErrParams = deepcopy(allSCIErrParams[(mpasFCRes,biasCorrectType.get(osName,None))])
475 
476  if osName is None or osName not in SCIErrParams:
477  _logger.error('osName not available in SCIErrParams => '+osName)
478  os._exit(1)
479 
480  # varName, ch = vu.splitIntSuffix(insituParameters[vu.selfDepValue])
481  varName, ch = vu.splitIntSuffix(insituParameters[vu.selfHofXValue])
482  STD0 = SCIErrParams[osName][(int(ch), SCISTDName)]['ERR'][0]
483  STD1 = SCIErrParams[osName][(int(ch), SCISTDName)]['ERR'][1]
484  SCI0 = SCIErrParams[osName][(int(ch), SCISTDName)]['X'][0]
485  SCI1 = SCIErrParams[osName][(int(ch), SCISTDName)]['X'][1]
486  slope = (STD1 - STD0) / (SCI1 - SCI0)
487 
488  belowramp = lessEqualBound(SCI, SCI0, False)
489  aboveramp = greatEqualBound(SCI, SCI1, False)
490  onramp = betweenBounds(SCI, SCI0, SCI1, False)
491 
492  BTerr = np.full_like(SCI, np.NaN)
493  BTerr[belowramp] = STD0
494  BTerr[onramp] = STD0 + slope * (SCI[onramp] - SCI0)
495  BTerr[aboveramp] = STD1
496 
497  # BTdep = dbVals[insituParameters[vu.selfDepValue]]
498  BTobs = dbVals[insituParameters[vu.selfObsValue]]
499  BTbak = dbVals[insituParameters[vu.selfHofXValue]]
500  BTdep = np.subtract(BTbak, BTobs)
501 
502  return np.divide(BTdep, BTerr)
503 
504 
506  def __init__(self):
507  super().__init__()
508  self.SCISCI = SCIOkamoto()
509  self.baseVarsbaseVarsbaseVars = pu.uniqueMembers(self.baseVarsbaseVarsbaseVars + self.SCISCI.baseVars)
510 
511  def evaluate(self, dbVals, insituParameters):
512  SCI = self.SCISCI.evaluate(dbVals, insituParameters)
513  return super().evaluate(dbVals, insituParameters, OkamotoMethod, SCI)
514 
515 
517  def __init__(self):
518  super().__init__()
519  self.SCISCI = ScaledSCIOkamoto()
520  self.baseVarsbaseVarsbaseVars = pu.uniqueMembers(self.baseVarsbaseVarsbaseVars + self.SCISCI.baseVars)
521 
522  def evaluate(self, dbVals, insituParameters):
523  SCI = self.SCISCI.evaluate(dbVals, insituParameters)
524  return super().evaluate(dbVals, insituParameters, ScaleOkamotoMethod, SCI)
525 
526 
528  def __init__(self):
529  super().__init__()
530  self.SCISCI = SCIModHarnisch()
531  self.baseVarsbaseVarsbaseVars = pu.uniqueMembers(self.baseVarsbaseVarsbaseVars + self.SCISCI.baseVars)
532 
533  def evaluate(self, dbVals, insituParameters):
534  SCI = self.SCISCI.evaluate(dbVals, insituParameters)
535  return super().evaluate(dbVals, insituParameters, ModHarnischMethod, SCI)
536 
537 
539  def __init__(self):
540  super().__init__()
542  self.baseVarsbaseVarsbaseVars = pu.uniqueMembers(self.baseVarsbaseVarsbaseVars + self.SCISCI.baseVars)
543 
544  def evaluate(self, dbVals, insituParameters):
545  SCI = self.SCISCI.evaluate(dbVals, insituParameters)
546  return super().evaluate(dbVals, insituParameters, ScaleModHarnischMethod, SCI)
547 
548 #TODO: use shapefiles/polygons to describe geographic regions instead of lat/lon boxes, e.g.,
549 #def outsideRegion(dbVals, REGION_NAME):
550 # Note: depending on shape file definitions, LON may need to be -180 to 180 instead of 0 to 360
551 #
552 # shp = READ_SHAPEFILE(REGION_NAME)
553 # lons = dbVals['longitdue']
554 # lats = dbVals['latitude']
555 # nlocs = len(lons)
556 # REGIONS = isinsideregions(lons, lats, shp)
557 # return REGIONS
558 
559 
560 #=========================================
561 # generic wrappers for ObsFunction classes
562 #=========================================
564  def __init__(self, baseVars):
565  self.baseVarsbaseVars = deepcopy(baseVars)
566  pass
567 
568  def dbVars(self, varName, fileFormat, outerIters_):
569  dbVars = []
570 
571  if (not isinstance(outerIters_, Iterable)
572  or isinstance(outerIters_,str)):
573  outerIters = [outerIters_]
574  else:
575  outerIters = outerIters_
576 
577  for baseVar in self.baseVarsbaseVars:
578  for outerIter in outerIters:
579  dbVar = vu.base2dbVar(
580  baseVar, varName, fileFormat, outerIter)
581  dbVars.append(dbVar)
582  return pu.uniqueMembers(dbVars)
583 
584 
586  def __init__(self, variable):
587  super().__init__([variable])
588 
589  def evaluate(self, dbVals, insituParameters):
590  return dbVals[insituParameters[self.baseVarsbaseVars[0]]]
591 
592 
594  def __init__(self, function):
595  self.functionfunction = function()
596  assert hasattr(self.functionfunction, 'baseVars'), \
597  ("ERROR, function class must have the baseVars attribute:", function)
598  super().__init__(self.functionfunction.baseVars)
599 
600  def evaluate(self, dbVals, insituParameters):
601  return self.functionfunction.evaluate(dbVals, insituParameters)
602 
603 
605  def __init__(self, config):
606  self.osNameosName = config['osName']
607  self.fileFormatfileFormat = config['fileFormat']
608 
609  variable = config['variable']
610  varIsString = isinstance(variable,str)
611  varIsClass = inspect.isclass(variable)
612  assert varIsString ^ varIsClass, \
613  ("ERROR: 'variable' must either be a String or a Class", config)
614 
615  if varIsString:
616  self.functionfunction = IdObsFunction(variable)
617 
618  if varIsClass:
619  self.functionfunction = ObsFunction(variable)
620 
621  def dbVars(self, varName, outerIters):
622  return self.functionfunction.dbVars(varName, self.fileFormatfileFormat, outerIters)
623 
624  def evaluate(self, dbVals, varName, outerIter):
625  # setup context-specific insituParameters for the evaluation
626  insituParameters = {}
627  for baseVar in self.functionfunction.baseVars:
628  insituParameters[baseVar] = vu.base2dbVar(
629  baseVar, varName, self.fileFormatfileFormat, outerIter)
630  insituParameters['osName'] = self.osNameosName
631 
632  # evaluate the function
633  self.resultresult = self.functionfunction.evaluate(dbVals, insituParameters)
634 
635 
636 #========================
637 # generic binning classes
638 #========================
639 
640 class BinFilter:
641  def __init__(self, config):
642  self.wherewhere = config['where']
643  tmp = config['bounds']
644  nBins = config['nBins']
645 
646  # note: for where functions that take multiple bounds (e.g., notEqualAnyBound)
647  # config['bounds'] can wrap an inner Iterable, e.g.,
648  # * config['bounds'] = [[0,1]] would apply the bounds 0 and 1 to all nBins bins
649  # * config['bounds'] = [[0,1], [5,6]] would apply the bounds 0 and 1 to the first bin,
650  # the bounds 5 and 6 to the second bin,
651  # and nBins must be equal to 2
652 
653  ibins = list(range(nBins))
654 
655  # allow for scalar and Iterable bounds
656  if (not isinstance(tmp, Iterable) or
657  isinstance(tmp, str)):
658  self.boundsbounds = np.empty(nBins, dtype=type(tmp))
659  for ii in ibins:
660  self.boundsbounds[ii] = tmp
661  else:
662  self.boundsbounds = np.empty(nBins, dtype=type(tmp[0]))
663 
664  # single element Iterable is applied uniformly to all bins
665  if len(tmp) == 1:
666  for ii in ibins:
667  self.boundsbounds[ii] = tmp[0]
668  # multiple element Iterable must be same length as nBins
669  elif len(tmp) == nBins:
670  for ii in ibins:
671  self.boundsbounds[ii] = tmp[ii]
672  else:
673  _logger.error("'bounds' must be a scalar, single-member Iterable, or an Iterable with the same length as 'values'!")
674  os._exit(1)
675  self.functionfunction = ObsFunctionWrapper(config)
676  self.except_diagsexcept_diags = config.get('except_diags', [])
677  self.mask_valuemask_value = config.get('mask_value', np.NaN)
678  #TODO: add other actions besides mask_value/exclude
679 
680 # def baseVars(self):
681 # return pu.uniqueMembers(self.function.baseVars)
682 
683  def dbVars(self, varName, outerIters):
684  dbVars = self.functionfunction.dbVars(
685  varName, outerIters)
686  return pu.uniqueMembers(dbVars)
687 
688  def evaluate(self, dbVals, varName, outerIter):
689  self.functionfunction.evaluate(dbVals, varName, outerIter)
690 
691  def apply(self, array, diagName, ibin):
692  newArray = deepcopy(array)
693 
694  if diagName not in self.except_diagsexcept_diags:
695  # remove locations where the mask is True
696  mask = self.wherewhere(self.functionfunction.result,(self.boundsbounds)[ibin])
697 
698  if len(mask) == len(newArray):
699  newArray[mask] = self.mask_valuemask_value
700  else:
701  _logger.error('BinFilter mask is incorrectly defined!')
702  os._exit(1)
703 
704  return newArray
705 
706 
707 exclusiveDiags = ['obs','bak','ana','SCI']
708 
709 class BinMethod:
710  def __init__(self, config):
711  #allows for scalar, str, and Iterable 'values'
712  tmp = config['values']
713  self.valuesvalues = []
714  if (not isinstance(tmp, Iterable) or
715  isinstance(tmp, str)):
716  self.valuesvalues += [tmp]
717  else:
718  self.valuesvalues += tmp
719 
720  self.excludeDiagsexcludeDiags = deepcopy(exclusiveDiags)
721  override = config.get('override_exclusiveDiags',[])
722  for diag in override:
723  if diag in self.excludeDiagsexcludeDiags:
724  self.excludeDiagsexcludeDiags.remove(diag)
725 
726  fconf = {}
727  fconf['osName'] = config['osName']
728  fconf['fileFormat'] = config['fileFormat']
729  fconf['nBins'] = len(self.valuesvalues)
730 
731  self.filtersfilters = []
732  for filterConf in config['filters']:
733  filterConf.update(fconf)
734  self.filtersfilters.append(BinFilter(filterConf))
735 
736  enoughBounds = False
737  for Filter in self.filtersfilters:
738  if len(Filter.bounds) == len(self.valuesvalues):
739  enoughBounds = True
740  assert enoughBounds, '\n\nERROR: BinMethod : at least one filter must have len(bounds) == len(values)!'
741 
742 # def baseVars(self):
743 # baseVars = []
744 # for Filter in self.filters:
745 # for variable in Filter.baseVars():
746 # baseVars.append(variable)
747 # return pu.uniqueMembers(baseVars)
748 
749  def dbVars(self, varName, outerIters=None):
750  dbVars = []
751  for Filter in self.filtersfilters:
752  dbVars += Filter.dbVars(
753  varName, outerIters)
754  return pu.uniqueMembers(dbVars)
755 
756  def evaluate(self, dbVals, varName, outerIter):
757  for ii in list(range(len(self.filtersfilters))):
758  self.filtersfilters[ii].evaluate(
759  dbVals, varName, outerIter)
760 
761  def apply(self, array, diagName, binVal):
762  ibin = self.valuesvalues.index(binVal)
763  masked_array = deepcopy(array)
764  for Filter in self.filtersfilters:
765  masked_array = Filter.apply(
766  masked_array, diagName, ibin)
767  return masked_array
np.maximum(zeros, np.subtract(BTclr, BTbak)))
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def dbVars(self, varName, fileFormat, outerIters_)
def __init__(self, baseVars)
def __init__(self, config)
def evaluate(self, dbVals, varName, outerIter)
def dbVars(self, varName, outerIters)
def apply(self, array, diagName, ibin)
def __init__(self, config)
def evaluate(self, dbVals, varName, outerIter)
def apply(self, array, diagName, binVal)
def dbVars(self, varName, outerIters=None)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def __init__(self, variable)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def __init__(self, function)
def dbVars(self, varName, outerIters)
def evaluate(self, dbVals, varName, outerIter)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters, SCISTDName, SCI)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def evaluate(self, dbVals, insituParameters)
def notEqualAnyBound(x, bounds, missingValue=True)
def lessEqualBound(x, bound, missingValue=True)
def greatEqualBound(x, bound, missingValue=True)
def betweenBounds(x, bound1, bound2, missingValue=True)
def equalBound(x, bound, missingValue=True)
def greatBound(x, bound, missingValue=True)
def lessBound(x, bound, missingValue=True)
def notEqualBound(x, bound, missingValue=True)