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

Last change on this file since 3696 was 3684, checked in by idelkadi, 4 years ago

Repository under svn of a first version of Multiatlas diagnostics for LMDZ. This version is adapted to be able to run a LMDZ multiatlas on an individual account on the ciclad machine of the IPSL. In this version, the parts to be modified are identified so as to subsequently adapt it to other machines.
This version is still under development.

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