source: LMDZ5/trunk/libf/cosp/phys_cosp.F90 @ 1925

Last change on this file since 1925 was 1925, checked in by idelkadi, 10 years ago
  • Re-ecriture des sorties du simulateur Cosp : Les 3 fichiers mensuel, journalier et haute frequence ne sont plus geres de facon classique par IOIPSL. Les fichiers "includes" ini_hist*COSP.h et write_hist*COSP.h sont supprimes et remplaces par 2 module :
    1. cosp_output_mod.F90 : ou sont crees les fichiers et ou sont definis les dimensions et les differents axes.
    2. cosp_output_write_mod.F90 : ou sont definis les variables diagnostiques a stocker dans ces fichiers et ou est geree leur ecriture.

L'utilisation d'XIO est prevu (A tester)

  • Rajouts d'autre variables diagnostiques : sunlit (=0 si nuit et =1 si jour) parasol_crefl (reflectance integree)
  • Correction du bug dans l'interface avec la physique : On ne veut pas que la distinction entre les nuages convectifs et stratiformes soient prise en compte dans les calculs Cosp. Dans l'interface avec la physique, sont passees en entree pour Cosp, les quantites totales (stratiforme + convective) du contenus en eau et d'autres variables. La fraction nuageuse convective calculee dans la physique est passee en entree pour Cosp dans la version buggee. Elle est remise a 0 dans cette version corrigee.
  • Mise a jour pour ISCCP : la fraction d'ensoleillement calculee par LMDZ est passee en argument pour Isccp
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 17.1 KB
Line 
1! Simulateur COSP : Cfmip Observation Simulator Package
2
3! ISCCP, Radar (QuickBeam), Lidar et Parasol (ACTSIM), MISR, RTTOVS
4!Idelkadi Abderrahmane Aout-Septembre 2009
5
6  subroutine phys_cosp( itap,dtime,freq_cosp, &
7                        ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
8                        ecrit_mth,ecrit_day,ecrit_hf, &
9                        Nptslmdz,Nlevlmdz,lon,lat, presnivs,overlaplmdz,sunlit, &
10                        ref_liq,ref_ice,fracTerLic,u_wind,v_wind,phis,phi,ph,p,skt,t, &
11                        sh,rh,tca,cca,mr_lsliq,mr_lsice,fl_lsrainI,fl_lssnowI, &
12                        fl_ccrainI,fl_ccsnowI,mr_ozone,dtau_s,dem_s)
13
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15!!!! Inputs :
16! itap,                                 !Increment de la physiq
17! dtime,                                !Pas de temps physiq
18! overlap,                              !Overlap type in SCOPS
19! Npoints,                              !Nb de points de la grille physiq
20! Nlevels,                              !Nb de niveaux verticaux
21! Ncolumns,                             !Number of subcolumns
22! lon,lat,                              !Longitudes et latitudes de la grille LMDZ
23! ref_liq,ref_ice,                      !Rayons effectifs des particules liq et ice (en microm)
24! fracTerLic,                           !Fraction terre a convertir en masque
25! u_wind,v_wind,                        !Vents a 10m ???
26! phi,                                  !Geopotentiel
27! phis,                                 !Geopotentiel sol
28! ph,                                   !pression pour chaque inter-couche
29! p,                                    !Pression aux milieux des couches
30! skt,t,                                !Temp au sol et temp 3D
31! sh,                                   !Humidite specifique
32! rh,                                   !Humidite relatif
33! tca,                                  !Fraction nuageuse
34! cca                                   !Fraction nuageuse convective
35! mr_lsliq,                             !Liq Cloud water content
36! mr_lsice,                             !Ice Cloud water content
37! mr_ccliq,                             !Convective Cloud Liquid water content 
38! mr_ccice,                             !Cloud ice water content
39! fl_lsrain,                            !Large scale precipitation lic
40! fl_lssnow,                            !Large scale precipitation ice
41! fl_ccrain,                            !Convective precipitation lic
42! fl_ccsnow,                            !Convective precipitation ice
43! mr_ozone,                             !Concentration ozone (Kg/Kg)
44! dem_s                                 !Cloud optical emissivity
45! dtau_s                                !Cloud optical thickness
46! emsfc_lw = 1.                         !Surface emissivity dans radlwsw.F90
47
48!!! Outputs :
49! calipso2D,                            !Lidar Low/heigh/Mean/Total-level Cloud Fraction
50! calipso3D,                            !Lidar Cloud Fraction (532 nm)
51! cfadlidar,                            !Lidar Scattering Ratio CFAD (532 nm)
52! parasolrefl,                          !PARASOL-like mono-directional reflectance
53! atb,                                  !Lidar Attenuated Total Backscatter (532 nm)
54! betamol,                              !Lidar Molecular Backscatter (532 nm)
55! cfaddbze,                             !Radar Reflectivity Factor CFAD (94 GHz)
56! clcalipso2,                           !Cloud frequency of occurrence as seen by CALIPSO but not CloudSat
57! dbze,                                 !Efective_reflectivity_factor
58! cltlidarradar,                        !Lidar and Radar Total Cloud Fraction
59! clMISR,                               !Cloud Fraction as Calculated by the MISR Simulator
60! clisccp2,                             !Cloud Fraction as Calculated by the ISCCP Simulator
61! boxtauisccp,                          !Optical Depth in Each Column as Calculated by the ISCCP Simulator
62! boxptopisccp,                         !Cloud Top Pressure in Each Column as Calculated by the ISCCP Simulator
63! tclisccp,                             !Total Cloud Fraction as Calculated by the ISCCP Simulator
64! ctpisccp,                             !Mean Cloud Top Pressure as Calculated by the ISCCP Simulator
65! tauisccp,                             !Mean Optical Depth as Calculated by the ISCCP Simulator
66! albisccp,                             !Mean Cloud Albedo as Calculated by the ISCCP Simulator
67! meantbisccp,                          !Mean all-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
68! meantbclrisccp                        !Mean clear-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
69!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
70
71  USE MOD_COSP_CONSTANTS
72  USE MOD_COSP_TYPES
73  USE MOD_COSP
74  USE mod_phys_lmdz_para
75  USE mod_grid_phy_lmdz
76  use ioipsl
77  use iophy
78  use cosp_output_mod
79  use cosp_output_write_mod
80 
81  IMPLICIT NONE
82
83  ! Local variables
84  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
85  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
86  character(len=512), save :: finput ! Input file name
87  character(len=512), save :: cmor_nl
88  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
89  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
90  integer, save :: Npoints      ! Number of gridpoints
91!$OMP THREADPRIVATE(Npoints)
92  integer, save :: Nlevels      ! Number of levels
93  Integer :: Nptslmdz,Nlevlmdz ! Nb de points issus de physiq.F
94  integer, save :: Nlr          ! Number of levels in statistical outputs
95  integer, save :: Npoints_it   ! Max number of gridpoints to be processed in one iteration
96  integer :: i
97  type(cosp_config),save :: cfg   ! Configuration options
98!$OMP THREADPRIVATE(cfg)
99  type(cosp_gridbox) :: gbx ! Gridbox information. Input for COSP
100  type(cosp_subgrid) :: sgx     ! Subgrid outputs
101  type(cosp_sgradar) :: sgradar ! Output from radar simulator
102  type(cosp_sglidar) :: sglidar ! Output from lidar simulator
103  type(cosp_isccp)   :: isccp   ! Output from ISCCP simulator
104  type(cosp_misr)    :: misr    ! Output from MISR simulator
105  type(cosp_vgrid)   :: vgrid   ! Information on vertical grid of stats
106  type(cosp_radarstats) :: stradar ! Summary statistics from radar simulator
107  type(cosp_lidarstats) :: stlidar ! Summary statistics from lidar simulator
108
109  integer :: t0,t1,count_rate,count_max
110  integer :: Nlon,Nlat,geomode
111  real,save :: radar_freq,k2,ZenAng,co2,ch4,n2o,co,emsfc_lw
112!$OMP THREADPRIVATE(emsfc_lw)
113  integer,dimension(RTTOV_MAX_CHANNELS),save :: Channels
114  real,dimension(RTTOV_MAX_CHANNELS),save :: Surfem
115  integer, save :: surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay
116  integer, save :: Nprmts_max_hydro,Naero,Nprmts_max_aero,lidar_ice_type
117  integer, save :: platform,satellite,Instrument,Nchannels
118  logical, save :: use_vgrid,csat_vgrid,use_precipitation_fluxes,use_reff
119
120! Declaration necessaires pour les sorties IOIPSL
121  integer :: ii
122  real    :: ecrit_day,ecrit_hf,ecrit_mth
123  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP
124
125  logical, save :: debut_cosp=.true.
126!$OMP THREADPRIVATE(debut_cosp)
127
128  include "dimensions.h"
129 
130!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
131  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
132  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
133                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
134                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
135  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
136  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit         
137  real,dimension(Nlevlmdz)        :: presnivs
138  integer                         :: itap,k,ip
139  real                            :: dtime,freq_cosp
140 
141!
142   namelist/COSP_INPUT/cmor_nl,overlap,isccp_topheight,isccp_topheight_direction, &
143              npoints_it,ncolumns,nlevels,use_vgrid,nlr,csat_vgrid,finput, &
144              radar_freq,surface_radar,use_mie_tables, &
145              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
146              lidar_ice_type,use_precipitation_fluxes,use_reff, &
147              platform,satellite,Instrument,Nchannels, &
148              Channels,Surfem,ZenAng,co2,ch4,n2o,co
149
150!---------------- End of declaration of variables --------------
151   
152
153!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
154! Read namelist with COSP inputs
155!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
156
157 if (debut_cosp) then
158  NPoints=Nptslmdz
159  Nlevels=Nlevlmdz
160 
161! Lecture du namelist input
162  CALL read_cosp_input
163
164! Clefs Outputs
165  call read_cosp_output_nl(cosp_output_nl,cfg)
166
167    if (overlaplmdz.ne.overlap) then
168       print*,'Attention overlaplmdz different de overlap lu dans namelist '
169    endif
170   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
171
172  print*,' Cles des differents simulateurs cosp :'
173  print*,' Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lrttov_sim', &
174          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lrttov_sim
175
176  endif ! debut_cosp
177
178!  print*,'Debut phys_cosp itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf ', &
179!          itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf
180!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
181! Allocate memory for gridbox type
182!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
183        print *, 'Allocating memory for gridbox type...'
184
185        call construct_cosp_gridbox(dble(itap),radar_freq,surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay,k2, &
186                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
187                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
188                                    use_precipitation_fluxes,use_reff, &
189                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
190                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
191       
192!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
193! Here code to populate input structure
194!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
195
196        print *, 'Populating input structure...'
197        gbx%longitude = lon
198        gbx%latitude = lat
199
200        gbx%p = p !
201        gbx%ph = ph
202        gbx%zlev = phi/9.81
203
204        zlev_half(:,1) = phis(:)/9.81
205        do k = 2, Nlevels
206          do ip = 1, Npoints
207           zlev_half(ip,k) = phi(ip,k)/9.81 + &
208               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
209          enddo
210        enddo
211        gbx%zlev_half = zlev_half
212
213        gbx%T = T
214        gbx%q = rh*100.
215        gbx%sh = sh
216! On ne veut pas que cosp distingue les nuages stratiformes et convectifs
217! on passe les contenus totaux (conv+strat)
218        gbx%cca = 0. !convective_cloud_amount (1)
219        gbx%tca = tca ! total_cloud_amount (1)
220        gbx%psfc = ph(:,1) !pression de surface
221        gbx%skt  = skt !Skin temperature (K)
222
223        do ip = 1, Npoints
224          if (fracTerLic(ip).ge.0.5) then
225             gbx%land(ip) = 1.
226          else
227             gbx%land(ip) = 0.
228          endif
229        enddo
230        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
231! A voir l equivalent LMDZ (u10m et v10m)
232        gbx%u_wind  = u_wind !eastward_wind (m s-1)
233        gbx%v_wind  = v_wind !northward_wind
234
235! sunlit calcule a partir de la fraction d ensoleillement par jour
236      do ip = 1, Npoints
237        if (sunlit(ip).le.0.) then
238           gbx%sunlit(ip)=0.
239        else
240           gbx%sunlit(ip)=1.
241        endif
242      enddo
243
244! A voir l equivalent LMDZ
245  mr_ccliq = 0.0
246  mr_ccice = 0.0
247        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
248        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
249        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
250        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
251! A revoir
252        fl_lsrain = fl_lsrainI + fl_ccrainI
253        fl_lssnow = fl_lssnowI + fl_ccsnowI
254        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
255        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
256!  A voir l equivalent LMDZ
257        fl_lsgrpl=0.
258        fl_ccsnow = 0.
259        fl_ccrain = 0.
260        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
261        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
262        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
263
264     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
265     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
266     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
267     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
268
269        ! ISCCP simulator
270        gbx%dtau_s   = dtau_s
271        gbx%dtau_c   = 0.
272        gbx%dem_s    = dem_s
273        gbx%dem_c    = 0.
274
275! Surafce emissivity
276       emsfc_lw = 1.
277               
278!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
279        ! Define new vertical grid
280!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
281        print *, 'Defining new vertical grid...'
282        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
283
284!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
285       ! Allocate memory for other types
286!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
287        print *, 'Allocating memory for other types...'
288        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
289        call construct_cosp_sgradar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
290        call construct_cosp_radarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
291        call construct_cosp_sglidar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
292        call construct_cosp_lidarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
293        call construct_cosp_isccp(cfg,Npoints,Ncolumns,Nlevels,isccp)
294        call construct_cosp_misr(cfg,Npoints,misr)
295
296!+++++++++++++ Open output files and define output files axis !+++++++++++++
297     if (debut_cosp) then
298
299        print *, ' Open outpts files and define axis'
300        call cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp, &
301                              ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, &
302                              ecrit_mth, ecrit_day, ecrit_hf, use_vgrid, vgrid)
303
304        debut_cosp=.false.
305      endif ! debut_cosp
306!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
307        ! Call simulator
308!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
309        print *, 'Calling simulator...'
310        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,stradar,stlidar)
311!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
312
313!!!!!!!!!!!!!!!!!! Ecreture des sorties Cosp !!!!!!!!!!!!!!r!!!!!!:!!!!!
314       print *, 'Calling write output'
315        call cosp_output_write(Nlevlmdz, Npoints, Ncolumns, itap, dtime, freq_COSP, cfg, gbx, sglidar, stlidar, isccp)
316
317!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
318        ! Deallocate memory in derived types
319!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
320        print *, 'Deallocating memory...'
321        call free_cosp_gridbox(gbx)
322        call free_cosp_subgrid(sgx)
323        call free_cosp_sgradar(sgradar)
324        call free_cosp_radarstats(stradar)
325        call free_cosp_sglidar(sglidar)
326        call free_cosp_lidarstats(stlidar)
327        call free_cosp_isccp(isccp)
328        call free_cosp_misr(misr)
329        call free_cosp_vgrid(vgrid) 
330 
331!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
332  ! Time in s. Only for testing purposes
333!  call system_clock(t1,count_rate,count_max)
334!  print *,(t1-t0)*1.0/count_rate
335 
336  CONTAINS
337 
338  SUBROUTINE read_cosp_input
339   
340    IF (is_master) THEN
341      OPEN(10,file=cosp_input_nl,status='old')
342      READ(10,nml=cosp_input)
343      CLOSE(10)
344    ENDIF
345    CALL bcast(cmor_nl)
346    CALL bcast(overlap)
347    CALL bcast(isccp_topheight)
348    CALL bcast(isccp_topheight_direction)
349    CALL bcast(npoints_it)
350    CALL bcast(ncolumns)
351    CALL bcast(nlevels)
352    CALL bcast(use_vgrid)
353    CALL bcast(nlr)
354    CALL bcast(csat_vgrid)
355    CALL bcast(finput)
356    CALL bcast(radar_freq)
357    CALL bcast(surface_radar)
358    CALL bcast(use_mie_tables)
359    CALL bcast(use_gas_abs)
360    CALL bcast(do_ray)
361    CALL bcast(melt_lay)
362    CALL bcast(k2)
363    CALL bcast(Nprmts_max_hydro)
364    CALL bcast(Naero)
365    CALL bcast(Nprmts_max_aero)
366    CALL bcast(lidar_ice_type)
367    CALL bcast(use_precipitation_fluxes)
368    CALL bcast(use_reff)
369    CALL bcast(platform)
370    CALL bcast(satellite)
371    CALL bcast(Instrument)
372    CALL bcast(Nchannels)
373    CALL bcast(Channels)
374    CALL bcast(Surfem)
375    CALL bcast(ZenAng)
376    CALL bcast(co2)
377    CALL bcast(ch4)
378    CALL bcast(n2o)
379    CALL bcast(co)
380!$OMP BARRIER 
381  END SUBROUTINE read_cosp_input
382
383end subroutine phys_cosp
Note: See TracBrowser for help on using the repository browser.