source: LMDZ6/trunk/libf/phylmd/cosp2/phys_cosp2.F90 @ 3364

Last change on this file since 3364 was 3364, checked in by idelkadi, 6 years ago

Implementation de COSPv2 dans LMDZ.

File size: 21.2 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 First Version
5!Idelkadi Abderrahmane Nov 2015 version v1.4.0
6
7  subroutine phys_cosp2( itap,dtime,freq_cosp, &
8                        ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
9                        ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
10                        Nptslmdz,Nlevlmdz,lon,lat, presnivs,overlaplmdz,sunlit, &
11                        ref_liq,ref_ice,fracTerLic,u_wind,v_wind,phis,phi,ph,p,skt,t, &
12                        sh,rh,tca,cca,mr_lsliq,mr_lsice,fl_lsrainI,fl_lssnowI, &
13                        fl_ccrainI,fl_ccsnowI,mr_ozone,dtau_s,dem_s)
14
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16!!!! Inputs :
17! itap,                                 !Increment de la physiq
18! dtime,                                !Pas de temps physiq
19! overlap,                              !Overlap type in SCOPS
20! Npoints,                              !Nb de points de la grille physiq
21! Nlevels,                              !Nb de niveaux verticaux
22! Ncolumns,                             !Number of subcolumns
23! lon,lat,                              !Longitudes et latitudes de la grille LMDZ
24! ref_liq,ref_ice,                      !Rayons effectifs des particules liq et ice (en microm)
25! fracTerLic,                           !Fraction terre a convertir en masque
26! u_wind,v_wind,                        !Vents a 10m ???
27! phi,                                  !Geopotentiel
28! phis,                                 !Geopotentiel sol
29! ph,                                   !pression pour chaque inter-couche
30! p,                                    !Pression aux milieux des couches
31! skt,t,                                !Temp au sol et temp 3D
32! sh,                                   !Humidite specifique
33! rh,                                   !Humidite relatif
34! tca,                                  !Fraction nuageuse
35! cca                                   !Fraction nuageuse convective
36! mr_lsliq,                             !Liq Cloud water content
37! mr_lsice,                             !Ice Cloud water content
38! mr_ccliq,                             !Convective Cloud Liquid water content 
39! mr_ccice,                             !Cloud ice water content
40! fl_lsrain,                            !Large scale precipitation lic
41! fl_lssnow,                            !Large scale precipitation ice
42! fl_ccrain,                            !Convective precipitation lic
43! fl_ccsnow,                            !Convective precipitation ice
44! mr_ozone,                             !Concentration ozone (Kg/Kg)
45! dem_s                                 !Cloud optical emissivity
46! dtau_s                                !Cloud optical thickness
47! emsfc_lw = 1.                         !Surface emissivity dans radlwsw.F90
48
49!!! Outputs :
50! calipso2D,                            !Lidar Low/heigh/Mean/Total-level Cloud Fraction
51! calipso3D,                            !Lidar Cloud Fraction (532 nm)
52! cfadlidar,                            !Lidar Scattering Ratio CFAD (532 nm)
53! parasolrefl,                          !PARASOL-like mono-directional reflectance
54! atb,                                  !Lidar Attenuated Total Backscatter (532 nm)
55! betamol,                              !Lidar Molecular Backscatter (532 nm)
56! cfaddbze,                             !Radar Reflectivity Factor CFAD (94 GHz)
57! clcalipso2,                           !Cloud frequency of occurrence as seen by CALIPSO but not CloudSat
58! dbze,                                 !Efective_reflectivity_factor
59! cltlidarradar,                        !Lidar and Radar Total Cloud Fraction
60! clMISR,                               !Cloud Fraction as Calculated by the MISR Simulator
61! clisccp2,                             !Cloud Fraction as Calculated by the ISCCP Simulator
62! boxtauisccp,                          !Optical Depth in Each Column as Calculated by the ISCCP Simulator
63! boxptopisccp,                         !Cloud Top Pressure in Each Column as Calculated by the ISCCP Simulator
64! tclisccp,                             !Total Cloud Fraction as Calculated by the ISCCP Simulator
65! ctpisccp,                             !Mean Cloud Top Pressure as Calculated by the ISCCP Simulator
66! tauisccp,                             !Mean Optical Depth as Calculated by the ISCCP Simulator
67! albisccp,                             !Mean Cloud Albedo as Calculated by the ISCCP Simulator
68! meantbisccp,                          !Mean all-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
69! meantbclrisccp                        !Mean clear-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
70
71!!! AI rajouter les nouvelles sorties
72!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
73
74!! AI rajouter
75!#include "cosp_defs.h"
76!  USE MOD_COSP_CONSTANTS
77!  USE MOD_COSP_TYPES
78!  USE MOD_COSP
79  USE mod_phys_lmdz_para
80  USE mod_grid_phy_lmdz
81  use ioipsl
82  use iophy
83  use cosp_output_mod
84  use cosp_output_write_mod
85!  use MOD_COSP_Modis_Simulator, only : cosp_modis
86!  use mod_cosp_config, only : vgrid_zl,vgrid_zu,vgrid_z
87  USE MOD_COSP_INTERFACE_v1p4, ONLY: cosp                   => cosp_interface_v1p4,      &
88                                     cosp_gridbox,construct_cosp_vgrid,                  &
89                                     construct_cosp_gridbox,                             &
90                                     free_cosp_gridbox      => destroy_cosp_gridbox,     &
91                                     free_cosp_sgradar      => destroy_cosp_sgradar,     &
92                                     free_cosp_radarstats   => destroy_cosp_radarstats,  &
93                                     free_cosp_sglidar      => destroy_cosp_sglidar,     &
94                                     free_cosp_lidarstats   => destroy_cosp_lidarstats,  &
95                                     free_cosp_isccp        => destroy_cosp_isccp,       &
96                                     free_cosp_misr         => destroy_cosp_misr,        &
97                                     free_cosp_rttov        => destroy_cosp_rttov,       &
98                                     free_cosp_modis        => destroy_cosp_modis,       &
99                                     free_cosp_vgrid        => destroy_cosp_vgrid,       &
100                                     free_cosp_subgrid      => destroy_cosp_subgrid,     &
101                                     construct_cosp_subgrid,cosp_config,cosp_subgrid,    &
102                                     cosp_sglidar,cosp_lidarstats,                       &
103                                     construct_cosp_lidarstats,construct_cosp_sglidar,   &
104                                     cosp_isccp,construct_cosp_isccp,cosp_misr,          &
105                                     construct_cosp_misr,cosp_rttov,construct_cosp_rttov,&
106                                     cosp_sgradar,cosp_radarstats,                       &
107                                     construct_cosp_radarstats,construct_cosp_sgradar,   &
108                                     cosp_modis,construct_cosp_modis,                    &
109                                     cosp_vgrid,I_CVCLIQ,I_LSCLIQ,I_CVCICE,I_LSCICE,     &
110                                     I_LSRAIN,I_LSSNOW,I_LSGRPL,I_CVRAIN,I_CVSNOW
111
112  use cosp_read_otputkeys
113
114  IMPLICIT NONE
115
116  ! Local variables
117  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
118  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
119  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
120  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
121  integer, save :: Npoints      ! Number of gridpoints
122!$OMP THREADPRIVATE(Npoints)
123  integer, save :: Nlevels      ! Number of levels
124  Integer :: Nptslmdz,Nlevlmdz ! Nb de points issus de physiq.F
125  integer, save :: Nlvgrid          ! Number of levels in statistical outputs
126  integer, save :: Npoints_it   ! Max number of gridpoints to be processed in one iteration
127!  integer :: i
128  type(cosp_config),save :: cfg   ! Configuration options
129!$OMP THREADPRIVATE(cfg)
130  type(cosp_gridbox) :: gbx ! Gridbox information. Input for COSP
131  type(cosp_subgrid) :: sgx     ! Subgrid outputs
132  type(cosp_sgradar) :: sgradar ! Output from radar simulator
133  type(cosp_sglidar) :: sglidar ! Output from lidar simulator
134  type(cosp_isccp)   :: isccp   ! Output from ISCCP simulator
135!! AI rajout modis
136  type(cosp_modis)   :: modis   ! Output from MODIS simulator
137!!
138  type(cosp_misr)    :: misr    ! Output from MISR simulator
139!! AI rajout rttovs
140  type(cosp_rttov)   :: rttov   ! Output from RTTOV
141!!
142  type(cosp_vgrid)   :: vgrid   ! Information on vertical grid of stats
143  type(cosp_radarstats) :: stradar ! Summary statistics from radar simulator
144  type(cosp_lidarstats) :: stlidar ! Summary statistics from lidar simulator
145
146  integer :: t0,t1,count_rate,count_max
147  integer :: Nlon,Nlat
148  real(wp),save :: radar_freq,k2,ZenAng,co2,ch4,n2o,co,emsfc_lw
149 
150!$OMP THREADPRIVATE(emsfc_lw)
151  integer,dimension(RTTOV_MAX_CHANNELS),save :: Channels
152  real(wp),dimension(RTTOV_MAX_CHANNELS),save :: Surfem
153  integer, save :: surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay
154  integer, save :: Nprmts_max_hydro,Naero,Nprmts_max_aero,lidar_ice_type
155  integer, save :: platform,satellite,Instrument,Nchannels
156  logical, save :: use_vgrid,csat_vgrid,use_precipitation_fluxes,use_reff
157
158! Declaration necessaires pour les sorties IOIPSL
159  integer :: ii
160  real    :: ecrit_day,ecrit_hf,ecrit_mth, missing_val
161  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, ok_all_xml
162
163  logical, save :: debut_cosp=.true.
164!$OMP THREADPRIVATE(debut_cosp)
165
166  logical, save :: first_write=.true.
167!$OMP THREADPRIVATE(first_write)
168
169!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
170  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
171  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
172                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
173                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
174  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
175  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit         
176  real,dimension(Nlevlmdz)        :: presnivs
177  integer                         :: itap,k,ip
178  real                            :: dtime,freq_cosp
179  real,dimension(2)               :: time_bnds
180
181  double precision                            :: d_dtime
182  double precision,dimension(2)               :: d_time_bnds
183
184   namelist/COSP_INPUT/overlap,isccp_topheight,isccp_topheight_direction, &
185              npoints_it,ncolumns,use_vgrid,Nlvgrid,csat_vgrid, &
186              radar_freq,surface_radar,use_mie_tables, &
187              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
188              lidar_ice_type,use_precipitation_fluxes,use_reff, &
189              platform,satellite,Instrument,Nchannels, &
190              Channels,Surfem,ZenAng,co2,ch4,n2o,co
191
192!---------------- End of declaration of variables --------------
193
194!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
195! Read namelist with COSP inputs
196!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
197
198 print*,'Entree phys_cosp2'
199 if (debut_cosp) then
200  NPoints=Nptslmdz
201  Nlevels=Nlevlmdz
202 
203! Lecture du namelist input
204!  CALL read_cosp_input
205   IF (is_master) THEN
206      OPEN(10,file=cosp_input_nl,status='old')
207      READ(10,nml=cosp_input)
208      CLOSE(10)
209   ENDIF
210!$OMP BARRIER
211    CALL bcast(overlap)
212    CALL bcast(isccp_topheight)
213    CALL bcast(isccp_topheight_direction)
214    CALL bcast(npoints_it)
215    CALL bcast(ncolumns)
216    CALL bcast(use_vgrid)
217    CALL bcast(Nlvgrid)
218    CALL bcast(csat_vgrid)
219    CALL bcast(radar_freq)
220    CALL bcast(surface_radar)
221    CALL bcast(use_mie_tables)
222    CALL bcast(use_gas_abs)
223    CALL bcast(do_ray)
224    CALL bcast(melt_lay)
225    CALL bcast(k2)
226    CALL bcast(Nprmts_max_hydro)
227    CALL bcast(Naero)
228    CALL bcast(Nprmts_max_aero)
229    CALL bcast(lidar_ice_type)
230    CALL bcast(use_precipitation_fluxes)
231    CALL bcast(use_reff)
232    CALL bcast(platform)
233    CALL bcast(satellite)
234    CALL bcast(Instrument)
235    CALL bcast(Nchannels)
236    CALL bcast(Channels)
237    CALL bcast(Surfem)
238    CALL bcast(ZenAng)
239    CALL bcast(co2)
240    CALL bcast(ch4)
241    CALL bcast(n2o)
242    CALL bcast(co)
243
244    print*,'ok read  cosp_input_nl'
245
246! Clefs Outputs initialisation
247  call cosp_outputkeys_init(cfg)
248!!!   call cosp_outputkeys_test(cfg)
249  print*,' Cles des differents simulateurs cosp a itap :',itap
250  print*,'Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lmodis_sim,Lrttov_sim,Lstats', &
251          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lmodis_sim, &
252          cfg%Lrttov_sim,cfg%Lstats
253
254    if (overlaplmdz.ne.overlap) then
255       print*,'Attention overlaplmdz different de overlap lu dans namelist '
256    endif
257   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
258
259  endif ! debut_cosp
260
261!!! Ici on modifie les cles logiques pour les outputs selon les champs actives dans les .xml
262  if ((itap.gt.1).and.(first_write))then
263#ifdef CPP_XIOS
264    call read_xiosfieldactive(cfg)
265#else
266    call read_cosp_output_nl(itap,cosp_output_nl,cfg)
267#endif
268    first_write=.false.
269
270    print*,' Cles des differents simulateurs cosp a itap :',itap
271    print*,'Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lmodis_sim,Lrttov_sim,Lstats', &
272          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lmodis_sim, &
273          cfg%Lrttov_sim,cfg%Lstats
274  endif
275
276  time_bnds(1) = dtime-dtime/2.
277  time_bnds(2) = dtime+dtime/2.
278
279  d_time_bnds=time_bnds
280  d_dtime=dtime
281
282!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
283! Allocate memory for gridbox type
284!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
285! AI mars 2017
286!        print *, 'Allocating memory for gridbox type...'
287
288! Surafce emissivity
289        emsfc_lw = 1.
290
291        call construct_cosp_gridbox(d_dtime,d_time_bnds,radar_freq,surface_radar,use_mie_tables,use_gas_abs, &
292                                    do_ray,melt_lay,k2, &
293                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
294                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
295                                    use_precipitation_fluxes,use_reff, &
296                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
297                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
298       
299!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
300! Here code to populate input structure
301!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
302
303!        print *, 'Populating input structure...'
304        gbx%longitude = lon
305        gbx%latitude = lat
306
307        gbx%p = p !
308        gbx%ph = ph
309        gbx%zlev = phi/9.81
310
311        zlev_half(:,1) = phis(:)/9.81
312        do k = 2, Nlevels
313          do ip = 1, Npoints
314           zlev_half(ip,k) = phi(ip,k)/9.81 + &
315               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
316          enddo
317        enddo
318        gbx%zlev_half = zlev_half
319
320        gbx%T = T
321        gbx%q = rh*100.
322        gbx%sh = sh
323! On ne veut pas que cosp distingue les nuages stratiformes et convectifs
324! on passe les contenus totaux (conv+strat)
325        gbx%cca = 0. !convective_cloud_amount (1)
326        gbx%tca = tca ! total_cloud_amount (1)
327        gbx%psfc = ph(:,1) !pression de surface
328        gbx%skt  = skt !Skin temperature (K)
329
330        do ip = 1, Npoints
331          if (fracTerLic(ip).ge.0.5) then
332             gbx%land(ip) = 1.
333          else
334             gbx%land(ip) = 0.
335          endif
336        enddo
337        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
338! A voir l equivalent LMDZ (u10m et v10m)
339        gbx%u_wind  = u_wind !eastward_wind (m s-1)
340        gbx%v_wind  = v_wind !northward_wind
341
342! sunlit calcule a partir de la fraction d ensoleillement par jour
343!      do ip = 1, Npoints
344!        if (sunlit(ip).le.0.) then
345!           gbx%sunlit(ip)=0.
346!        else
347!           gbx%sunlit(ip)=1.
348!        endif
349!      enddo
350       gbx%sunlit=sunlit
351
352! A voir l equivalent LMDZ
353  mr_ccliq = 0.0
354  mr_ccice = 0.0
355        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
356        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
357        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
358        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
359! A revoir
360        fl_lsrain = fl_lsrainI + fl_ccrainI
361        fl_lssnow = fl_lssnowI + fl_ccsnowI
362        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
363        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
364!  A voir l equivalent LMDZ
365        fl_lsgrpl=0.
366        fl_ccsnow = 0.
367        fl_ccrain = 0.
368        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
369        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
370        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
371
372     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
373     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
374!! AI A revoir
375     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
376     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
377
378        ! ISCCP simulator
379        gbx%dtau_s   = dtau_s
380        gbx%dtau_c   = 0.
381        gbx%dem_s    = dem_s
382        gbx%dem_c   = 0.
383
384!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
385        ! Define new vertical grid
386!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
387!        print *, 'Defining new vertical grid...'
388        call construct_cosp_vgrid(gbx,Nlvgrid,use_vgrid,csat_vgrid,vgrid)
389        print*,'Nlvgrid,use_vgrid,csat_vgrid ',Nlvgrid,use_vgrid,csat_vgrid
390        print*,'vgrid%z, vgrid%mz ',vgrid%z, vgrid%mz
391!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
392       ! Allocate memory for other types
393!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
394!        print *, 'Allocating memory for other types...'
395        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
396        if (cfg%Lradar_sim) call construct_cosp_sgradar(Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
397        if (cfg%Lradar_sim) call construct_cosp_radarstats(Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
398        if (cfg%Llidar_sim) call construct_cosp_sglidar(Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
399        if (cfg%Llidar_sim) call construct_cosp_lidarstats(Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
400        if (cfg%Lisccp_sim) call construct_cosp_isccp(Npoints,Ncolumns,Nlevels,isccp)
401!! AI rajout
402        if (cfg%Lmodis_sim) call construct_cosp_modis(Npoints,modis)
403!!
404        if (cfg%Lmisr_sim) call construct_cosp_misr(Npoints,misr)
405!        call construct_cosp_rttov(cfg,Npoints,Nchannels,rttov)
406
407    if (debut_cosp) then
408      !$OMP MASTER
409        print *, ' Open outpts files and define axis'
410        call cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp, &
411                              ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml, &
412                              ecrit_mth, ecrit_day, ecrit_hf, vgrid)
413      !$OMP END MASTER
414      !$OMP BARRIER
415       debut_cosp=.false.
416      endif ! debut_cosp
417
418!!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
419        ! Call simulator
420!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
421      print*,'call simulateur'
422       call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,   &
423                                 isccp,misr,modis,rttov,stradar,stlidar)
424
425
426!!!!!!!!!!!!!!!!!! Ecreture des sorties Cosp !!!!!!!!!!!!!!r!!!!!!:!!!!!
427
428       print *, 'Calling write output'
429        call cosp_output_write(Nlevlmdz, Npoints, Ncolumns, itap, dtime, freq_COSP, missing_val, &
430                               cfg, gbx, vgrid, sglidar, sgradar, stlidar, stradar, &
431                               isccp, misr, modis)
432!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
433        ! Deallocate memory in derived types
434!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
435      print *, 'Deallocating memory...'
436     call free_cosp_gridbox(gbx)
437      print *, 'free_cosp_gridbox'
438     call free_cosp_subgrid(sgx)
439      print *, 'ffree_cosp_subgrid'
440     call free_cosp_vgrid(vgrid)
441      print *, 'free_cosp_vgrid'
442     if (cfg%Lradar_sim) call free_cosp_sgradar(sgradar)
443     if (cfg%Lradar_sim) call free_cosp_radarstats(stradar)
444     if (cfg%Llidar_sim) call free_cosp_sglidar(sglidar)
445     if (cfg%Llidar_sim) call free_cosp_lidarstats(stlidar)
446     if (cfg%Lisccp_sim) call free_cosp_isccp(isccp)
447     if (cfg%Lmisr_sim)  call free_cosp_misr(misr)
448     if (cfg%Lmodis_sim) call free_cosp_modis(modis)
449     if (cfg%Lrttov_sim) call free_cosp_rttov(rttov)
450     print *, 'End phys_cosp2'
451!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
452  ! Time in s. Only for testing purposes
453!  call system_clock(t1,count_rate,count_max)
454!  print *,(t1-t0)*1.0/count_rate
455 
456end subroutine phys_cosp2
Note: See TracBrowser for help on using the repository browser.