source: BOL/Multi_atlas/atlas/lmdz_SE.py.ref @ 5272

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