MPAS-JEDI
GenerateABEIFactors.py
Go to the documentation of this file.
1 #!/usr/bin/env python3
2 
3 import GenerateABEIFactorsArgs
4 
5 from basic_plot_functions import plotDistri
6 import binning_utils as bu
7 import predefined_configs as pconf
8 from collections import defaultdict
9 from collections.abc import Iterable
10 import config as conf
11 from copy import deepcopy
12 import diag_utils as du
13 import fnmatch
14 import glob
15 from JediDB import JediDB
16 from JediDBArgs import obsFKey
17 import Interpolate
18 import logging
19 import logsetup
20 import modelsp_utils as modelUtils
21 import multiprocessing as mp
22 from netCDF4 import Dataset
23 import numpy as np
24 import os
25 #from scipy.spatial import cKDTree
26 import stat_utils as su
27 import var_utils as vu
28 
29 _logger = logging.getLogger(__name__)
30 
32  '''
33  Diagnose observation-space adaptive background error inflation
34  Driven by
35  - static selections in conf
36  - command-line arguments in GenerateABEIFactorsArgs
37  '''
38  def __init__(self):
39  global argsProcessor
40  self.namename = 'GenerateABEIFactors'
41  self.argsargs = GenerateABEIFactorsArgs.args
42  self.loggerlogger = logging.getLogger(self.namename)
43 
44  # localization length scale for inflation projection
45  self.localizationLengthlocalizationLength = 2.0e5 #meters ... same as all-sky gaussian thinning distance
46  #self.localizationLength = 7.2e5 #meters ... 6x grid spacing
47  #self.localizationLength = 2.e6 #meters same as background ensemble localization length
48 
49  # construct mean DB into 0th member slot
50  self.loggerlogger.info('database path: '+self.argsargs.dbPath)
51  self.dbdb = JediDB(self.argsargs.dbPath)
52  self.osNamesosNames = self.argsargs.IRInstruments.split(',')
53  dbOSNames = list(self.dbdb.ObsSpaceName.values())
54  self.osKeysosKeys = []
55  for inst in self.osNamesosNames:
56  assert inst in dbOSNames, 'GenerateABEIFactors::__init__ '+inst+' not included in JediDB'
57  for osKey, osName in self.dbdb.ObsSpaceName.items():
58  if inst == osName:
59  self.osKeysosKeys.append(osKey)
60  break
61 
62  def diagnose(self):
63  '''
64  conducts diagnoseObsSpace across multiple ObsSpaces in parallel
65  '''
66  # Loop over all experiment+observation combinations (keys) alphabetically
67  for kk, osKey in enumerate(self.osKeysosKeys):
68  self.loggerlogger.info(osKey)
69  #Note: processing must be done serially due to the nature of the inflation
70  # update procedure... could process in parallel across channels/obsVars
71  self.diagnoseObsSpacediagnoseObsSpacediagnoseObsSpace(self.dbdb, osKey, {'resetLambda': kk==0})
72 
73  def getObsVars(self, db, osKey):
74  # create observed variable list by selecting those variables in the
75  # obs feedback files (obsFKey) with the proper suffix
76 # if self.args.jediAppName == 'variational':
77 # markerSuffix = vu.depbgGroup
78 # elif self.args.jediAppName == 'hofx':
79 # markerSuffix = vu.hofxGroup
80 # else:
81 # logger.error('JEDI Application is not supported:: '+self.args.jediAppName)
82 #
83 # return db.varList(osKey, obsFKey, markerSuffix)
84  return [
85  'brightness_temperature_8',
86  'brightness_temperature_9',
87  'brightness_temperature_10',
88  ]
89  @staticmethod
91  binVarConfigs = defaultdict(list)
92  binVarConfigs[vu.obsVarQC] += [bu.goodQCMethod]
93  return binVarConfigs
94 
95  @staticmethod
96  def getBinMethods(binVarConfigs, ObsSpaceName, fileFormat):
97  binMethods = {}
98  for binVarKey, binMethodKeys in binVarConfigs.items():
99  binVarConfig = pconf.binVarConfigs.get(binVarKey,pconf.nullBinVarConfig)
100  for binMethodKey in binMethodKeys:
101  config = binVarConfig.get(binMethodKey,pconf.nullBinMethod).copy()
102 
103  if (len(config['values']) < 1 or
104  len(config['filters']) < 1): continue
105 
106  config['osName'] = ObsSpaceName
107  config['fileFormat'] = fileFormat
108  binMethods[(binVarKey,binMethodKey)] = bu.BinMethod(config)
109 
110  return binMethods
111 
112  @staticmethod
114  return ['ABEILambda']
115 
116  @staticmethod
117  def getDBVars(obsVars, binMethods, diagnosticConfigs):
118  dbVars = []
119  for varName in obsVars:
120  for diagName, diagnosticConfig in diagnosticConfigs.items():
121  if 'ObsFunction' not in diagnosticConfig: continue
122 
123  # variables for diagnostics
124  for grpVar in diagnosticConfig['ObsFunction'].dbVars(
125  varName, diagnosticConfig['outerIter']):
126  if diagnosticConfig[vu.mean]:
127  dbVars.append(grpVar)
128 
129  # variables for binning
130  for (binVarKey,binMethodKey), binMethod in binMethods.items():
131  for grpVar in binMethod.dbVars(
132  varName, diagnosticConfig['outerIter']):
133  dbVars.append(grpVar)
134 
135  dbVars.append(vu.latMeta)
136  dbVars.append(vu.lonMeta)
137 
138  return dbVars
139 
140 
141  def diagnoseObsSpace(self, db, osKey, localConfig):
142  # osKey - key of db dict to process
143  logger = logging.getLogger(self.namename+'.diagnoseObsSpace('+osKey+')')
144  osName = self.osNamesosNames[self.osKeysosKeys.index(osKey)]
145 
146  # initialize mean db file handles
147  db.initHandles(osKey)
148 
149  ###############################################
150  ## Extract constructor info about the ObsSpace
151  ###############################################
152 
153  ObsSpaceName = db.ObsSpaceName[osKey]
154  ObsSpaceInfo = conf.DiagSpaceConfig[ObsSpaceName]
155  ObsSpaceGrp = ObsSpaceInfo['DiagSpaceGrp']
156 
157  obsVars = self.getObsVarsgetObsVarsgetObsVars(db, osKey)
158 
159  ########################################################
160  ## Construct dictionary of binMethods for this ObsSpace
161  ########################################################
162 
163  binVarConfigs = self.getBinVarConfigsgetBinVarConfigsgetBinVarConfigs()
164  binMethods = self.getBinMethodsgetBinMethodsgetBinMethods(binVarConfigs, ObsSpaceName, db.fileFormat(osKey, obsFKey))
165 
166  ######################################
167  ## Construct diagnostic configurations
168  ######################################
169 
170  diagnosticConfigs = du.diagnosticConfigs(
171  self.getDiagNamesgetDiagNamesgetDiagNames(), ObsSpaceName,
172  includeEnsembleDiagnostics = False,
173  fileFormat = db.fileFormat(osKey, obsFKey))
174 
175 # TODO: change above line after bugfix/hofxDiagsCleanup is merged into develop
176 # self.getDiagNames(), ObsSpaceName, False)
177 
178  ######################################
179  ## Generate list of required variables
180  ######################################
181 
182  dbVars = self.getDBVarsgetDBVarsgetDBVars(obsVars, binMethods, diagnosticConfigs)
183 
184 
185  ##################################
186  ## Read required variables from db
187  ##################################
188 
189  # read mean database variable values into memory
190  dbVals = db.readVars(osKey, dbVars)
191 
192  # destroy mean file handles
193  db.destroyHandles(osKey)
194 
195 
196  #################################################
197  ## Calculate and plot diagnostics for all obsVars
198  ## Calculate and write out inflation factors
199  #################################################
200 
201  # read model variables used for inflation calculation/templating
202  (modelLatsDeg, modelLonsDeg, EarthSphereR) = modelUtils.readGrid(gridFile=self.argsargs.modelGridFile, returnR=True)
203  modelLats = np.multiply(modelLatsDeg, np.pi / 180.0)
204  modelLons = np.multiply(modelLonsDeg, np.pi / 180.0)
205 
206  modelTemplateVars = {}
207  for templateType, templateVar in modelUtils.templateVariables.items():
208  modelTemplateVars[templateType] = {
209  'values': modelUtils.varRead(templateVar, self.argsargs.modelGridFile),
210  'dims': modelUtils.varDims(templateVar, self.argsargs.modelGridFile),
211  }
212 
213  # setup observation lat/lon
214  obsLatsDeg = dbVals[vu.latMeta]
215  obsLonsDeg = dbVals[vu.lonMeta]
216  obsLats = np.multiply(obsLatsDeg, np.pi / 180.0)
217  obsLons = np.multiply(obsLonsDeg, np.pi / 180.0)
218  obsnLocs = len(obsLons)
219 
220  # setup interpolation object
221  model2obs = Interpolate.InterpolateLonLat(
222  modelLons, modelLats,
223  weightMethod = 'barycentric',
224  Radius = EarthSphereR)
225  model2obs.initWeights(obsLons, obsLats)
226 
227  #TODO: move these attributes to diag_utils
228  diagColors = {}
229  minmaxValue = {}
230  diagColors['ABEILambda'] = 'BuPu'
231  minmaxValue['ABEILambda'] = [
232  bu.ABEILambda().minLambda,
233  bu.ABEILambda().maxLambda,
234  ]
235 
236  for diagName, diagnosticConfig in diagnosticConfigs.items():
237  if 'ObsFunction' not in diagnosticConfig: continue
238 
239  logger.info('Calculating/writing diagnostic stats for:')
240  logger.info('Diagnostic => '+diagName)
241  Diagnostic = diagnosticConfig['ObsFunction']
242  outerIter = diagnosticConfig['outerIter']
243 
244  for varName in obsVars:
245  logger.info('Variable => '+varName)
246 
247  varShort, varUnits = vu.varAttributes(varName)
248 
249  Diagnostic.evaluate(dbVals, varName, outerIter)
250  diagValues = Diagnostic.result
251  nLocs = len(diagValues)-np.isnan(diagValues).sum()
252 
253  if nLocs == 0:
254  logger.warning('All missing values for diagnostic: '+diagName)
255 
256  for (binVarKey,binMethodKey), binMethod in binMethods.items():
257  if diagName in binMethod.excludeDiags: continue
258 
259  binVarName, binGrpName = vu.splitObsVarGrp(binVarKey)
260  binVarShort, binVarUnits = vu.varAttributes(binVarName)
261 
262  # initialize binMethod filter function result
263  # NOTE: binning is performed using mean values
264  # and not ensemble member values
265  binMethod.evaluate(dbVals, varName, outerIter)
266 
267  for binVal in binMethod.values:
268  # apply binMethod filters for binVal
269  binnedDiagnostic = binMethod.apply(diagValues,diagName,binVal)
270 
271  # Plot horizontal distribution of binnedDiagnostic
272  if self.argsargs.plotLambda:
273  logger.info('plotting obs-space diagnostic')
274  dotsize = 9.0
275  color = diagColors[diagName]
276  minmax = minmaxValue.get(diagName, [None, None])
277  nLocs = len(binnedDiagnostic)-np.isnan(binnedDiagnostic).sum()
278  plotDistri(
279  obsLatsDeg, obsLonsDeg, binnedDiagnostic,
280  osName, varShort, varUnits, osName+'_'+binVarName+'='+binVal,
281  nLocs, diagName,
282  minmax[0], minmax[1], dotsize, color)
283 
284  #TODO: create a base class to be used by both this and DiagnoseObsStatistics
285  #TODO: encapsulate below behavior in a member function
286  if diagName == 'ABEILambda':
287  lambdaVarName = 'modelLambda_'+varShort
288  resetLambda = localConfig.get('resetLambda', True)
289  inflationFile = varShort+'_'+self.argsargs.inflationFile
290 
291  if not os.path.exists(inflationFile) or resetLambda:
292  # create new NC file with same header info as grid file
293  logger.info('creating fresh inflation file')
294  modelUtils.createHeaderOnlyFile(
295  self.argsargs.modelGridFile,
296  inflationFile,
297  date = self.argsargs.datetime,
298  )
299 
300  if modelUtils.hasVar(lambdaVarName, inflationFile):
301  # read previous values for lambdaVarName (i.e., from previous sensor)
302  modelLambda = modelUtils.varRead(lambdaVarName, inflationFile)
303  else:
304  modelLambda = np.full_like(modelTemplateVars['1D-c']['values'],
305  bu.ABEILambda().minLambda)
306 
307  ## Project observation-space inflation factors to
308  # model-space using localization function
309  # in sequential loop over the observation locations
310  # Minamide and Zhang (2018), eq. 10
311  logger.info('projecting inflation factors to model space')
312  for obsInd, (obsLambda, obsLat, obsLon) in enumerate(list(zip(
313  binnedDiagnostic, obsLats, obsLons))):
314  if np.isfinite(obsLambda):
315  gcDistances = Interpolate.HaversineDistance(
316  obsLon, obsLat,
317  modelLons, modelLats,
318  EarthSphereR)
319  scales = gcDistances / self.localizationLengthlocalizationLength
320  rho = np.zeros(scales.shape)
321 
322  #See Greybush et al. (2011) for discussion
323  # https://doi.org/10.1175/2010MWR3328.1
324  # Greybush attributes Gaspari and Cohn (1999) with prescribing
325  # the Gaussian localization function in model-space and
326  # relates it to Hunt et al. (2007) application to obs-space
327  # ... could revisit rho function for more optimal results
328  crit = scales < 3.5
329  rho[crit] = np.exp(- (scales[crit] ** 2) / 2.0)
330  crit = (np.abs(rho) > 0.0)
331  modelLocInds = np.arange(0, len(rho))[crit]
332 
333  # interpolate current lambda field to this obs location
334  interpLambda = model2obs.applyAtIndex(modelLambda, obsInd)
335 
336  # update lambda field
337  modelLambda[modelLocInds] = modelLambda[modelLocInds] + \
338  rho[crit] * (obsLambda - interpLambda)
339 
340  # truncate negative valuess caused by obs-space discretization
341  modelLambda[modelLambda < bu.ABEILambda().minLambda] = \
342  bu.ABEILambda().minLambda
343 
344  if self.argsargs.plotLambda:
345  # Plot horizontal distribution of modelLambda
346  logger.info('plotting model-space inflation factors')
347  dotsize = 9.0
348  color = diagColors[diagName]
349  minmax = minmaxValue.get(diagName, [None, None])
350  plotDistri(
351  modelLatsDeg, modelLonsDeg, modelLambda,
352  'model-'+osName.upper()+' '+diagName, varShort, '',
353  'model-'+osName+'_'+binVarName+'='+binVal,
354  0, diagName,
355  minmax[0], minmax[1], dotsize, color)
356 
357  logger.info('writing inflation factors to NC file')
358 
359  # write universal inflation factor used for all state variables below
360  attrs = {
361  'units': 'unitless',
362  'long_name': 'Column-wise Inflation factor derived from '+varShort,
363 
364  }
365  modelUtils.varWrite(lambdaVarName, modelLambda, inflationFile,
366  attrs, modelTemplateVars['1D-c']['dims'], 'double')
367 
368  # write state-variable-specific inflation factors
369  for modelInflateVarName in modelUtils.inflationVariables:
370  templateInfo = modelUtils.variableTraits[modelInflateVarName]
371  if modelUtils.hasVar(modelInflateVarName, self.argsargs.modelGridFile):
372  templateVarName = modelInflateVarName
373  templateVar = {
374  'attrs': modelUtils.varAttrs(modelInflateVarName, self.argsargs.modelGridFile),
375  }
376  else:
377  templateVarName = modelUtils.templateVariables[templateInfo['templateVar']]
378  templateVar = {
379  'attrs': {
380  'units': templateInfo['units'],
381  'long_name': templateInfo['long_name'],
382  },
383  }
384  templateVar['values'] = modelUtils.varRead(templateVarName, self.argsargs.modelGridFile)
385  templateVar['dims'] = modelUtils.varDims(templateVarName, self.argsargs.modelGridFile)
386  templateVar['datatype'] = modelUtils.varDatatype(templateVarName, self.argsargs.modelGridFile)
387 
388  shape = templateVar['values'].shape
389  nDims = len(shape)
390  if nDims == 1:
391  modelInflateVar = modelLambda
392  elif nDims == 2:
393  modelInflateVar = np.empty_like(templateVar['values'])
394  for level in np.arange(0, shape[1]):
395  modelInflateVar[:,level] = modelLambda
396 
397  modelUtils.varWrite(modelInflateVarName, modelInflateVar, inflationFile,
398  templateVar['attrs'], templateVar['dims'], templateVar['datatype'],
399  )
400 
401  #END binMethod.values LOOP
402  #END binMethods tuple LOOP
403  #END obsVars LOOP
404  #END diagnosticConfigs LOOP
405 
406 
407 #=========================================================================
408 # main program
409 def main():
410  _logger.info('Starting '+__name__)
411 
412  abei = GenerateABEIFactors()
413 
414  # diagnose abei
415  abei.diagnose()
416 
417  _logger.info('Finished '+__name__+' successfully')
418 
419 if __name__ == '__main__': main()
def getBinMethods(binVarConfigs, ObsSpaceName, fileFormat)
def diagnoseObsSpace(self, db, osKey, localConfig)
def getDBVars(obsVars, binMethods, diagnosticConfigs)
def HaversineDistance(lon1, lat1, lon2, lat2, R=1.)
Definition: Interpolate.py:770
def plotDistri(lats, lons, values, ObsType, VarName, var_unit, out_name, nstation, levbin, dmin=None, dmax=None, dotsize=6, color="rainbow")