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

Last change on this file since 2137 was 2137, checked in by idelkadi, 10 years ago

Implementation de XIOS pour les sorties du simulateur COSP

  • 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, ok_all_xml, &
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  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
87  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
88  integer, save :: Npoints      ! Number of gridpoints
89!$OMP THREADPRIVATE(Npoints)
90  integer, save :: Nlevels      ! Number of levels
91  Integer :: Nptslmdz,Nlevlmdz ! Nb de points issus de physiq.F
92  integer, save :: Nlr          ! Number of levels in statistical outputs
93  integer, save :: Npoints_it   ! Max number of gridpoints to be processed in one iteration
94  integer :: i
95  type(cosp_config),save :: cfg   ! Configuration options
96!$OMP THREADPRIVATE(cfg)
97  type(cosp_gridbox) :: gbx ! Gridbox information. Input for COSP
98  type(cosp_subgrid) :: sgx     ! Subgrid outputs
99  type(cosp_sgradar) :: sgradar ! Output from radar simulator
100  type(cosp_sglidar) :: sglidar ! Output from lidar simulator
101  type(cosp_isccp)   :: isccp   ! Output from ISCCP simulator
102  type(cosp_misr)    :: misr    ! Output from MISR simulator
103  type(cosp_vgrid)   :: vgrid   ! Information on vertical grid of stats
104  type(cosp_radarstats) :: stradar ! Summary statistics from radar simulator
105  type(cosp_lidarstats) :: stlidar ! Summary statistics from lidar simulator
106
107  integer :: t0,t1,count_rate,count_max
108  integer :: Nlon,Nlat,geomode
109  real,save :: radar_freq,k2,ZenAng,co2,ch4,n2o,co,emsfc_lw
110!$OMP THREADPRIVATE(emsfc_lw)
111  integer,dimension(RTTOV_MAX_CHANNELS),save :: Channels
112  real,dimension(RTTOV_MAX_CHANNELS),save :: Surfem
113  integer, save :: surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay
114  integer, save :: Nprmts_max_hydro,Naero,Nprmts_max_aero,lidar_ice_type
115  integer, save :: platform,satellite,Instrument,Nchannels
116  logical, save :: use_vgrid,csat_vgrid,use_precipitation_fluxes,use_reff
117
118! Declaration necessaires pour les sorties IOIPSL
119  integer :: ii
120  real    :: ecrit_day,ecrit_hf,ecrit_mth
121  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, ok_all_xml
122
123  logical, save :: debut_cosp=.true.
124!$OMP THREADPRIVATE(debut_cosp)
125
126  include "dimensions.h"
127 
128!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
129  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
130  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
131                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
132                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
133  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
134  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit         
135  real,dimension(Nlevlmdz)        :: presnivs
136  integer                         :: itap,k,ip
137  real                            :: dtime,freq_cosp
138 
139!
140   namelist/COSP_INPUT/overlap,isccp_topheight,isccp_topheight_direction, &
141              npoints_it,ncolumns,use_vgrid,nlr,csat_vgrid, &
142              radar_freq,surface_radar,use_mie_tables, &
143              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
144              lidar_ice_type,use_precipitation_fluxes,use_reff, &
145              platform,satellite,Instrument,Nchannels, &
146              Channels,Surfem,ZenAng,co2,ch4,n2o,co
147
148!---------------- End of declaration of variables --------------
149   
150
151!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
152! Read namelist with COSP inputs
153!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
154
155 if (debut_cosp) then
156  NPoints=Nptslmdz
157  Nlevels=Nlevlmdz
158 
159! Lecture du namelist input
160  CALL read_cosp_input
161
162! Clefs Outputs
163  call read_cosp_output_nl(cosp_output_nl,cfg)
164
165    if (overlaplmdz.ne.overlap) then
166       print*,'Attention overlaplmdz different de overlap lu dans namelist '
167    endif
168   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
169
170  print*,' Cles des differents simulateurs cosp :'
171  print*,' Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lrttov_sim', &
172          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lrttov_sim
173
174  endif ! debut_cosp
175
176!  print*,'Debut phys_cosp itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf ', &
177!          itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf
178!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179! Allocate memory for gridbox type
180!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
181        print *, 'Allocating memory for gridbox type...'
182
183        call construct_cosp_gridbox(dble(itap),radar_freq,surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay,k2, &
184                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
185                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
186                                    use_precipitation_fluxes,use_reff, &
187                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
188                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
189       
190!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
191! Here code to populate input structure
192!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
193
194        print *, 'Populating input structure...'
195        gbx%longitude = lon
196        gbx%latitude = lat
197
198        gbx%p = p !
199        gbx%ph = ph
200        gbx%zlev = phi/9.81
201
202        zlev_half(:,1) = phis(:)/9.81
203        do k = 2, Nlevels
204          do ip = 1, Npoints
205           zlev_half(ip,k) = phi(ip,k)/9.81 + &
206               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
207          enddo
208        enddo
209        gbx%zlev_half = zlev_half
210
211        gbx%T = T
212        gbx%q = rh*100.
213        gbx%sh = sh
214! On ne veut pas que cosp distingue les nuages stratiformes et convectifs
215! on passe les contenus totaux (conv+strat)
216        gbx%cca = 0. !convective_cloud_amount (1)
217        gbx%tca = tca ! total_cloud_amount (1)
218        gbx%psfc = ph(:,1) !pression de surface
219        gbx%skt  = skt !Skin temperature (K)
220
221        do ip = 1, Npoints
222          if (fracTerLic(ip).ge.0.5) then
223             gbx%land(ip) = 1.
224          else
225             gbx%land(ip) = 0.
226          endif
227        enddo
228        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
229! A voir l equivalent LMDZ (u10m et v10m)
230        gbx%u_wind  = u_wind !eastward_wind (m s-1)
231        gbx%v_wind  = v_wind !northward_wind
232
233! sunlit calcule a partir de la fraction d ensoleillement par jour
234      do ip = 1, Npoints
235        if (sunlit(ip).le.0.) then
236           gbx%sunlit(ip)=0.
237        else
238           gbx%sunlit(ip)=1.
239        endif
240      enddo
241
242! A voir l equivalent LMDZ
243  mr_ccliq = 0.0
244  mr_ccice = 0.0
245        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
246        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
247        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
248        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
249! A revoir
250        fl_lsrain = fl_lsrainI + fl_ccrainI
251        fl_lssnow = fl_lssnowI + fl_ccsnowI
252        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
253        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
254!  A voir l equivalent LMDZ
255        fl_lsgrpl=0.
256        fl_ccsnow = 0.
257        fl_ccrain = 0.
258        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
259        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
260        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
261
262     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
263     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
264     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
265     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
266
267        ! ISCCP simulator
268        gbx%dtau_s   = dtau_s
269        gbx%dtau_c   = 0.
270        gbx%dem_s    = dem_s
271        gbx%dem_c    = 0.
272
273! Surafce emissivity
274       emsfc_lw = 1.
275               
276!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
277        ! Define new vertical grid
278!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
279        print *, 'Defining new vertical grid...'
280        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
281
282!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
283       ! Allocate memory for other types
284!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
285        print *, 'Allocating memory for other types...'
286        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
287        call construct_cosp_sgradar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
288        call construct_cosp_radarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
289        call construct_cosp_sglidar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
290        call construct_cosp_lidarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
291        call construct_cosp_isccp(cfg,Npoints,Ncolumns,Nlevels,isccp)
292        call construct_cosp_misr(cfg,Npoints,misr)
293
294!+++++++++++++ Open output files and define output files axis !+++++++++++++
295     if (debut_cosp) then
296
297      !$OMP MASTER
298        print *, ' Open outpts files and define axis'
299        call cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp, &
300                              ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml, &
301                              ecrit_mth, ecrit_day, ecrit_hf, use_vgrid, vgrid)
302      !$OMP END MASTER
303      !$OMP BARRIER
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, &
316                               cfg, gbx, vgrid, sglidar, stlidar, isccp)
317
318!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
319        ! Deallocate memory in derived types
320!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
321        print *, 'Deallocating memory...'
322        call free_cosp_gridbox(gbx)
323        call free_cosp_subgrid(sgx)
324        call free_cosp_sgradar(sgradar)
325        call free_cosp_radarstats(stradar)
326        call free_cosp_sglidar(sglidar)
327        call free_cosp_lidarstats(stlidar)
328        call free_cosp_isccp(isccp)
329        call free_cosp_misr(misr)
330        call free_cosp_vgrid(vgrid) 
331 
332!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
333  ! Time in s. Only for testing purposes
334!  call system_clock(t1,count_rate,count_max)
335!  print *,(t1-t0)*1.0/count_rate
336 
337  CONTAINS
338 
339  SUBROUTINE read_cosp_input
340   
341    IF (is_master) THEN
342      OPEN(10,file=cosp_input_nl,status='old')
343      READ(10,nml=cosp_input)
344      CLOSE(10)
345    ENDIF
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(use_vgrid)
352    CALL bcast(nlr)
353    CALL bcast(csat_vgrid)
354    CALL bcast(radar_freq)
355    CALL bcast(surface_radar)
356    CALL bcast(use_mie_tables)
357    CALL bcast(use_gas_abs)
358    CALL bcast(do_ray)
359    CALL bcast(melt_lay)
360    CALL bcast(k2)
361    CALL bcast(Nprmts_max_hydro)
362    CALL bcast(Naero)
363    CALL bcast(Nprmts_max_aero)
364    CALL bcast(lidar_ice_type)
365    CALL bcast(use_precipitation_fluxes)
366    CALL bcast(use_reff)
367    CALL bcast(platform)
368    CALL bcast(satellite)
369    CALL bcast(Instrument)
370    CALL bcast(Nchannels)
371    CALL bcast(Channels)
372    CALL bcast(Surfem)
373    CALL bcast(ZenAng)
374    CALL bcast(co2)
375    CALL bcast(ch4)
376    CALL bcast(n2o)
377    CALL bcast(co)
378!$OMP BARRIER 
379  END SUBROUTINE read_cosp_input
380
381end subroutine phys_cosp
Note: See TracBrowser for help on using the repository browser.