source: BOL/Multi_atlas/atlas/lmdz_SE.py @ 4639

Last change on this file since 4639 was 4330, checked in by musat, 2 years ago

Ajout environnement et corrections

  • champ eva(poration) moyennes zonales
  • pour les appels a des fonctions climaf via mcdo (timavg, sqrt)

IonelaMusat?

File size: 17.7 KB
Line 
1# -*- coding: iso-8859-1 -*-
2# Created : S.Sénési - nov 2015
3
4"""
5Base pour le traitement des données SE de LMDZ :
6  - variables, saisons et grilles habituellement gérées,
7  - organisation des données SE
8  - renommage et re-scaling des variables
9  - dictionnaire de fichiers d'exemples pour les grilles traitées
10et une fonction pour rendre un objet CliMAF pour tuple (simu, variable, saison, grille)
11"""
12
13from climaf.api import *
14
15# Liste des variables LMDZ objet d'une cmorisation ou d'une projection+moyenne saisonniere
16#########################################################################################
17"""
18variables_list=['sfcWind','hfls','hfss','huss','hurs','pr',\
19    'rldscs','rlds','rlus','rsdscs','rsds','rsuscs','rsus','rsutcs','rsut','rlut','rlutcs','rsdt',\
20    'sfcWind','tas','ts','tauu','tauv','psl','zg500','hfns','ta','ua','va','wap','hus','hur',\
21    'cllcalipso','clmcalipso','clhcalipso','cltcalipso','clt','prw','tasmax','tasmin']
22"""
23variables_list=['sfcWind','hfls','hfss','bil','huss','hurs','pr',\
24    'rldscs','rlds','rlus','rsdscs','rsds','rsuscs','rsus','rsutcs','rsut','rlut','rlutcs','rsdt',\
25    'sfcWind','tas','sst','tauu','tauv','psl','zg500','hfns','ta','ua','va','wap','hus','hur',\
26    'cllcalipso','clmcalipso','clhcalipso','cltcalipso','clt','prw','tasmax','tasmin']
27
28#Dictionnaire des grilles connues et de fichiers d'exemple
29#########################################################################################
30grid_remapfiles={ 'VLR' : '/data/ssenesi/stable/VLR.nc'}
31
32# Saisons gérées et leur traduction CDO
33#########################################################################################
34#seasons={ 'YEAR':"", 'DJF' : " selmon,1,2,12", 'JJA' : " selmon,6,7,8" , 'JJAS' : "-timavg -selmon,6,7,8,9" }
35#seasons={ 'YEAR':"", 'DJF' : "-timavg -selmon,1,2,12", 'JJA' : "-timavg -selmon,6,7,8" , 'JJAS' : "-timavg -selmon,6,7,8,9" }
36seasons={ 'YEAR':"", 'DJF' : "timavg -selmon,1,2,12", 'JJA' : "timavg -selmon,6,7,8" , 'JJAS' : "timavg -selmon,6,7,8,9", 'MAM' : "timavg -selmon,3,4,5", 'SON' : "timavg -selmon,9,10,11" }
37
38
39#########################################################################################
40# Definitions pour l'acces aux obs gérées par Ionela
41#########################################################################################
42cproject('OBS', ('period','fx'),   ('root','/thredds/ipsl/fabric/lmdz/SE/CMOR/OBS'))  #, '/data/musat/LMDZ6/SE/CMOR/OBS'))
43dataloc(project='OBS',url='${root}/${variable}.nc')
44# Définiton de variables dérivées pour ce projet
45
46#derive('OBS','rstt','minus','rsdt','rsut')
47#derive('OBS','rstscs','minus','rsdscs','rsuscs')
48#derive('OBS','rsts','minus','rsds','rsus')
49
50#derive('OBS','cress','minus','rsds','rsdscs')
51#derive('OBS','crels','minus','rlds','rldscs')
52#derive('OBS','crets','plus','crels','cress')
53#derive('OBS','crest','minus','rsutcs','rsut')
54#derive('OBS','crelt','minus','rlutcs','rlut')
55#derive('OBS','crett','plus','crest','crelt')
56
57# Variables dérivées
58derive('OBS','rstt','minus','rsdt','rsut')
59derive('OBS','rsts','minus','rsds','rsus')
60#BUH here: IMderive('OBS','rlts','minus','rlds','rlus')
61derive('OBS','rlts','minus','rlus','rlds')
62derive('OBS','rtt','minus','rstt','rlut')
63derive('OBS','rts','plus','rsts','rlts')
64
65derive('OBS','rstscs','minus','rsdscs','rsuscs')
66derive('OBS','rsttcs','minus','rsdt'  ,'rsutcs')
67
68derive('OBS','cress','minus','rsds','rsdscs')
69derive('OBS','crels','minus','rlds','rldscs')
70derive('OBS','crets','plus','cress','crels')
71
72derive('OBS','crest','minus','rsutcs','rsut')
73derive('OBS','crelt','minus','rlutcs','rlut')
74derive('OBS','crett','plus','crest','crelt')
75
76derive('OBS','hfns','plus','hfls','hfss')
77derive('OBS','bil' ,'minus','rts','hfns')
78derive('OBS','tsmtas','minus','ts','tas')
79
80derive('OBS','rlah','minus','rlts','rlut')
81derive('OBS','rtmp','plus','rldscs','rlutcs')
82derive('OBS','rlahcs','minus','rlus','rtmp')
83
84derive('OBS','rlahcre','minus','rlah','rlahcs')
85
86derive('OBS','rsah','minus','rstt','rsts')
87derive('OBS','rsahcs','minus','rsttcs','rstscs')
88derive('OBS','rsahcre','minus','rsah','rsahcs')
89
90derive('OBS','rah','plus','rsah','rlah')
91derive('OBS','rahcs','plus','rsahcs','rlahcs')
92derive('OBS','rahcre','minus','rah','rahcs')
93
94cscript('mask_val_inf_1','cdo setrtomiss,0,1 ${in} ${out}')
95derive('OBS','rsdt_mask','mask_val_inf_1','rsdt')
96derive('OBS','albt','divide','rsut','rsdt_mask')
97derive('OBS','albtcs','divide','rsutcs','rsdt_mask')
98
99
100derive('OBS','albs','divide','rsus','rsds')
101derive('OBS','dtr','minus','tasmax','tasmin')
102
103calias('OBS', 'hflsevap', 'hfls', scale=1./2.5e6)
104
105#########################################################################################
106# Définition logique des donnees de type 'SE' de l'IPSL, et de leur localisation
107# La periode n'est pas traitee comme un objet 'period' de CliMAF, mais comme une chaine
108# de caracteres dont le nom de facette est 'years'
109
110# Le nommage des données découvert chez F.Hourdin est traité. On utilise la facette 'root'
111# pour indiquer le répertoire de base
112# Les re-scaling habituels sont decrits a l'aide de la fonction calias
113# Les variables dérivées sont décrites aussi
114#########################################################################################
115
116cproject('SE', 
117         ('frequency','seasonnal'), 
118         ('period','fx'), 
119         'years',
120         ('root','/thredds/ipsl/musat/lmdz/SE/ORIG'), 
121         separator='|')
122#exemple de nom de fichier : NPv3.1ada_SE_1982_1991_1M_histmthCOSP.nc
123pattern='${root}/${simulation}_SE_${years}_1M_histmth*.nc'
124dataloc(project='SE',url=pattern, organization='generic')
125   
126   
127calias('SE','hfls','flat',scale=-1.)
128calias('SE','hfss','sens',scale=-1.)
129calias('SE','pr','precip')
130calias('SE','sfcWind','wind10m')
131calias('SE','rldscs','LWdnSFCclr')
132calias('SE','rlds','LWdnSFC')
133calias('SE','rlus','LWupSFC')
134calias('SE','rsdscs','SWdnSFCclr')
135calias('SE','rsds','SWdnSFC')
136calias('SE','rsuscs','SWupSFCclr')
137calias('SE','rsus','SWupSFC')
138calias('SE','rsutcs','SWupTOAclr')
139calias('SE','rsut','SWupTOA')
140calias('SE','rsdt','SWdnTOA')
141calias('SE','rlut','topl')
142calias('SE','rlutcs','topl0')
143calias('SE','sfcWind','wind10m')
144calias('SE','tas','t2m')
145calias('SE','tasmax','ave_t2m_daily_max')
146calias('SE','tasmin','ave_t2m_daily_min')
147calias('SE','ts','tsol')
148calias('SE','sst','tsol_oce')
149calias('SE','huss','q2m')
150calias('SE','hurs','rh2m')
151calias('SE','tauu','taux_oce')
152calias('SE','tauv','tauy_oce')
153#calias('SE','psl','slp')
154calias('SE','psl','slp', scale=0.01)
155calias('SE','zg500','z500')
156#calias('SE','pslhPa','slp', scale=0.01)
157#calias('OBS','pslhPa','psl', scale=0.01)
158
159#calias('SE','hfns','bils')
160calias('SE','ta','temp')
161calias('SE','ua','vitu')
162calias('SE','va','vitv')
163calias('SE','wap','vitw')
164calias('SE','hus','ovap')
165calias('SE','hur','rhum', scale=100.)
166calias('SE','clt','cldt', scale=100.)
167calias('SE','cltcalipso','cltcalipso', scale=100.)
168calias('SE','clhcalipso','clhcalipso', scale=100.)
169calias('SE','clmcalipso','clmcalipso', scale=100.)
170calias('SE','cllcalipso','cllcalipso', scale=100.)
171
172# Variables dérivées
173derive('SE','rstt','minus','rsdt','rsut')
174derive('SE','rsts','minus','rsds','rsus')
175#BUG HERE derive('SE','rlts','minus','rlds','rlus')
176derive('SE','rlts','minus','rlus','rlds')
177derive('SE','rltscs','minus','rldscs','rluscs')
178
179derive('SE','rtt','minus','rstt','rlut')
180derive('SE','rts','plus','rsts','rlts')
181
182derive('SE','rstscs','minus','rsdscs','rsuscs')
183derive('SE','rsttcs','minus','rsdt'  ,'rsutcs')
184
185derive('SE','cress','minus','rsds','rsdscs')
186derive('SE','crels','minus','rlds','rldscs')
187derive('SE','crets','plus','cress','crels')
188
189derive('SE','crest','minus','rsutcs','rsut')
190derive('SE','crelt','minus','rlutcs','rlut')
191derive('SE','crett','plus','crest','crelt')
192
193derive('SE','hfns','plus','hfls','hfss')
194derive('SE','bil' ,'minus','rts','hfns')
195derive('SE','tsmtas','minus','ts','tas')
196
197derive('SE','rlah','minus','rlts','rlut')
198derive('SE','rtmp','plus','rldscs','rlutcs')
199derive('SE','rlahcs','minus','rlus','rtmp')
200
201#let rlahcs=rlus-rldscs-rlutcs
202derive('SE','rlahcre','minus','rlah','rlahcs')
203#
204derive('SE','rsah','minus','rstt','rsts')
205derive('SE','rsahcs','minus','rsttcs','rstscs')
206derive('SE','rsahcre','minus','rsah','rsahcs')
207
208derive('SE','rah','plus','rsah','rlah')
209derive('SE','rahcs','plus','rsahcs','rlahcs')
210derive('SE','rahcre','minus','rah','rahcs')
211
212cscript('mask_valLMDZ_inf_1','cdo setrtomiss,0,1 ${in} ${out}')
213derive('SE','rsdt_maskLMDZ','mask_valLMDZ_inf_1','rsdt')
214derive('SE','albt','divide','rsut','rsdt_maskLMDZ')
215derive('SE','albtcs','divide','rsutcs','rsdt_maskLMDZ')
216#derive('SE','albt','divide','rsut','rsdt')
217#derive('SE','albtcs','divide','rsutcs','rsdt')
218
219derive('SE','albs','divide','rsus','rsds')
220
221derive('SE','dtr','minus','tasmax','tasmin')
222
223# -- P - E
224calias('SE', 'hflsevap', 'flat', scale=-1./2.5e6, filenameVar='histmth')
225derive('SE', 'pme', 'minus', 'pr' ,'hflsevap')
226
227
228# -- Atmospheric Variables on vertical levels
229for tmpvar in ['ua', 'va', 'ta', 'hus', 'hur']:
230    for tmplev in ['850', '700', '500', '200']:
231        derive('OBS', tmpvar+tmplev, 'ccdo', tmpvar, operator='intlevel,'+tmplev)
232        derive('SE', tmpvar+tmplev, 'ccdo', tmpvar, operator='intlevel,'+tmplev+'00')
233
234
235# let albt=100*(rsut/rsdt)
236# let albs=100*(rsus/rsds)
237# let tsk=ts+273.18
238# let pslhPa=psl/100.
239# let tasc=tas-273.16
240
241def all_SE_simulations():
242    """
243    Listage de toutes les simulations du projet par listage de leurs
244    fichiers de données et décodage de leur nom
245   
246    Il faudrait quelque chose de moins adhoc ...
247    """
248    import re
249    simus=[]
250    # Listage des fichiers de donnees
251    a=ds(project='SE',years="*", simulation='*', variable='*')
252    for f in a.baseFiles().split(' ') : 
253        basename=f.split('/')[-1].replace('_SE','')
254        #basename=re.sub(r'_1M_histmth.*.nc','',basename)
255        basename=re.sub(r'_1M_histmth*.*.nc','',basename)
256        if basename not in simus  : simus.append(basename)
257    return simus
258
259
260def svsg(simulation,variable,season='YEAR',grid='',root='/thredds/ipsl/musat/lmdz/SE/ORIG'):
261    """
262    Rend l'objet CliMAF pour une variable d'une simulation, pour les  données
263    - SE d'un intervalle d'annéees (forme YYY1_YYY2)
264    - ou OBS
265    et une saison et une grille
266
267    (Le nom de la fonction est la concaténation des intiales des arguments)
268
269    Exemples :
270    >>> svsg('NPv3.1ada','1982_1991','hurs')
271    >>> svsg(simulation='NPv3.1ada',years='1982_1991',variable='hurs',season='DJF, grid='')
272    >>> svsg('OBS','','hurs','JJA')
273
274    Il faut au préalable avoir déclaré les projets 'SE' et 'OBS'
275    """
276    if simulation != 'OBS' :
277        # Il faut identifier les annees dans le nom de la simu
278        yeardeb=simulation.split('_')[1]
279        yearfin=simulation.split('_')[2]
280        years=yeardeb+"_"+yearfin
281        simulation=simulation.split('_')[0]
282        data=ds(project='SE',variable=variable,years=years, simulation=simulation)
283    else:
284        data=ds(project='OBS',variable=variable)
285    if season != 'YEAR' : 
286        seas=ccdo(data,operator=seasons[season])
287    else : 
288        seas=time_average(data)
289    rds=seas
290    if grid != '' : 
291        if grid not in grid_remapfiles : 
292            print("La grille %s n'est pas connue"%options.grid) ; 
293            return None
294        rds=regrid(seas,fds(grid_remapfiles[grid],period='fx'))
295    return rds
296
297def bias_SE(simu,variable,season='YEAR',grid='',root='/thredds/ipsl/musat/lmdz/SE/ORIG') :
298    """
299    Calcule le biais pour une variable d'une simu, vs 'OBS'
300    Rend -999. en cas de probleme
301    """
302    try :
303        sim=svsg(simu ,variable,season,grid,root)
304        ref=svsg('OBS',variable,season,grid)
305        if (grid == '' ) : ref=regrid(ref,sim)
306        dif=minus(sim,ref)
307        return cvalue(space_average(dif))
308    except :
309        return -999.
310
311def is3d(variable) :
312    if variable in ['ta','ua','va','hus','hur','wap','cl','clw','cli','mc','tro3'] :
313        return True
314    return False
315
316
317def varlongname(variable,short=False):
318   """
319   Returns the long name of variable
320   """
321   shortvarname = longvarname = variable
322   if variable=='tas':
323        longvarname  = '2M Temperature'
324        shortvarname = '2M Temp.'
325   if variable=='tasmax':
326        longvarname  = '2M daymax Temperature'
327        shortvarname = '2M daymx Temp.'
328   if variable=='tasmin':
329        longvarname  = '2M daymin Temperature'
330        shortvarname = '2M daymn Temp.'
331   if variable=='DRT':
332        longvarname  = '2M diurnal amplitude temperature'
333        shortvarname = '2M dtr '
334   if variable=='pr':
335        longvarname  = 'Precipitation'
336        shortvarname = 'Precip.'
337   if variable=='psl':
338        longvarname  = 'Sea Level Pressure'
339        shortvarname = 'Sea Level Pres.'
340   if variable=='ua':
341        longvarname  = 'Zonal Wind'
342        shortvarname = 'U-Wind'
343   if variable=='va':
344        longvarname  = 'Meridional Wind'
345        shortvarname = 'V-Wind'
346   if variable=='ta':
347        longvarname  = 'Air Temperature'
348        shortvarname = 'Air Temp.'
349   if variable=='hus':
350        longvarname  = 'Specific Humidity'
351        shortvarname = 'Sp. Humidity'
352   if variable=='hur':
353        longvarname  = 'Relative Humidity'
354        shortvarname = 'Rel. Humidity'
355   if variable=='huss':
356        longvarname  = 'Specific Humidity at Surface'
357        shortvarname = 'Sp. Humidity (surf)'
358   if variable=='rstt':
359        longvarname  = 'Rad SW Total TOA'
360        shortvarname = 'Rad SW Total TOA'
361   if variable=='rsts':
362        longvarname  = 'Total SW rad surface'
363        shortvarname = 'Total SW rad surf.'
364   if variable=='rtt':
365        longvarname  = 'Total Radiation TOA'
366        shortvarname = 'Total Rad. TOA'
367   if variable=='crelt':
368        longvarname  = 'Longwave Cloud Radiative Effect TOA'
369        shortvarname = 'LW CRE TOA'
370   if variable=='crest':
371        longvarname  = 'Shortwave Cloud Radiative Effect TOA'
372        shortvarname = 'SW CRE TOA'
373   if variable=='crett':
374        longvarname  = 'Total CRE TOA'
375        shortvarname = 'Total CRE TOA'
376   if variable=='hfls':
377        longvarname  = 'Latent Heat Flux'
378        shortvarname = 'Latent HF'
379   if variable=='hfss':
380        longvarname  = 'Sensible Heat Flux'
381        shortvarname = 'Sensible HF'
382   if variable=='hfns':
383        longvarname  = 'Surface Total Heat Flux'
384        shortvarname = 'Surf. Tot. HF'
385   if variable=='zg500':
386        longvarname  = '500mb geopotential height'
387        shortvarname = ''
388   if variable=='rsut':
389        longvarname  = 'Upward SW rad TOA'
390        shortvarname = ''
391   if variable=='rlut':
392        longvarname  = 'Outgoing Long Wave Radiation'
393        shortvarname = 'OLR'
394   if variable=='rlutcs':
395        longvarname  = 'Clear Sky OLR'
396        shortvarname = 'Clear sky OLR'
397   if variable=='albs':
398        longvarname  = 'Surface albedo'
399        shortvarname = ''
400   if variable=='albt':
401        longvarname  = 'Planetary albedo'
402        shortvarname = ''
403   if variable=='albtcs':
404        longvarname  = 'Clear Sky Planetary albedo'
405        shortvarname = ''
406   if variable=='cress':
407        longvarname  = 'SW CRE surface'
408        shortvarname = ''
409   if variable=='crels':
410        longvarname  = 'LW CRE surface'
411        shortvarname = ''
412   if variable=='crets':
413        longvarname  = 'Total CRE surface'
414        shortvarname = ''
415   if variable=='rts':
416        longvarname  = 'Total radiation surface'
417        shortvarname = ''
418   if variable=='rah':
419        longvarname  = 'Atm. Rad. Heat.'
420        shortvarname = ''
421   if variable=='rahcs':
422        longvarname  = 'Atm. Rad. Heat. - clear sky'
423        shortvarname = ''
424   if variable=='rahcre':
425        longvarname  = 'Atm. rad. Heat. - CRE'
426        shortvarname = ''
427   if variable=='rsah':
428        longvarname  = 'Atm. SW Heat.'
429        shortvarname = ''
430   if variable=='rsahcs':
431        longvarname  = 'Atm. SW Heat. - Clear sky'
432        shortvarname = ''
433   if variable=='rsahcre':
434        longvarname  = 'Atm. SW Heat. - CRE'
435        shortvarname = ''
436   if variable=='rlah':
437        longvarname  = 'Atm. LW Heat.'
438        shortvarname = ''
439   if variable=='rlahcs':
440        longvarname  = 'Atm. LW Heat. - Clear sky'
441        shortvarname = ''
442   if variable=='rlahcre':
443        longvarname  = 'Atm. LW Heat. - CRE'
444        shortvarname = ''
445   if variable=='cltcalipso':
446        longvarname  = 'Total Cloud Cover'
447        shortvarname = ''
448   if variable=='cllcalipso':
449        longvarname  = 'Low Cloud Cover'
450        shortvarname = ''
451   if variable=='clmcalipso':
452        longvarname  = 'Medium Cloud Cover'
453        shortvarname = ''
454   if variable=='clhcalipso':
455        longvarname  = 'High Cloud Cover'
456        shortvarname = ''
457   if variable=='rlds':
458        longvarname  = 'Downward LW rad at Surface'
459        shortvarname = ''
460   if variable=='rldscs':
461        longvarname  = 'Downward LW rad at Surface - Clear Sky'
462        shortvarname = ''
463   if variable=='hurs':
464        longvarname  = 'Relative Humidity at Surface'
465        shortvarname = ''
466   if variable=='rlus':
467        longvarname  = 'Upward LW rad at Surface'
468        shortvarname = ''
469   if variable=='rsus':
470        longvarname  = 'Upward SW rad at Surface'
471        shortvarname = ''
472   if variable=='rsuscs':
473        longvarname  = 'Upward SW rad at Surface - Clear Sky'
474        shortvarname = ''
475   if variable=='rsdscs':
476        longvarname  = 'Downward SW rad at Surface - Clear Sky'
477        shortvarname = ''
478   if variable=='rsds':
479        longvarname  = 'Downward SW rad at Surface'
480        shortvarname = ''
481   if variable=='rsucs':
482        longvarname  = 'Upward SW rad at Surface - Clear Sky'
483        shortvarname = ''
484   if variable=='rsutcs':
485        longvarname  = 'Upward SW rad at TOA - Clear Sky'
486        shortvarname = ''
487   if variable=='tauu':
488        longvarname  = 'Zonal Wind Stress'
489        shortvarname = ''
490   if variable=='tauv':
491        longvarname  = 'Meridional Wind Stress'
492        shortvarname = ''
493   if variable=='pme':
494        longvarname  = 'P-E Precip-Evap(hfls/2.5e6) mm/day'
495        shortvarname = ''
496   #
497   if short:
498        return shortvarname
499   else:
500        return longvarname
501
502
503
Note: See TracBrowser for help on using the repository browser.