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

Last change on this file since 4228 was 3745, 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.
Update atlas scripts

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