source: BOL/Multi_atlas/Utils/lmdz_SE0.py @ 5452

Last change on this file since 5452 was 4335, checked in by musat, 2 years ago

Autre oubli
Ionela

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