source: LMDZ5/trunk/libf/phylmd/cosp/phys_cosp.F90 @ 2345

Last change on this file since 2345 was 2345, checked in by Ehouarn Millour, 9 years ago

Physics/dynamics separation:

  • move test_disvert_m to dynlonlat_phylonlat/phylmd since it is only used by ce0l and relies on dynamics.
  • put "config_inca" in tracinca_mod so physics routines can get the info from there rather than from control_mod.
  • get rid of references to "control_mod" from within the physics.

EM

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