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

Last change on this file since 4666 was 4619, checked in by yann meurdesoif, 12 months ago

Suppress usage of preprocessing key CPP_XIOS.
Wrapper file is used to suppress XIOS symbol when xios is not linked and not used (-io ioipsl)
The CPP_XIOS key is replaced in model by "using_xios" boolean variable to switch between IOIPSL or XIOS output.

YM

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 lmdz_xios, ONLY : using_xios
86!  use MOD_COSP_Modis_Simulator, only : cosp_modis
87!  use mod_cosp_config, only : vgrid_zl,vgrid_zu,vgrid_z
88  USE MOD_COSP_INTERFACE_v1p4, ONLY: cosp                   => cosp_interface_v1p4,      &
89                                     cosp_gridbox,construct_cosp_vgrid,                  &
90                                     construct_cosp_gridbox,                             &
91                                     free_cosp_gridbox      => destroy_cosp_gridbox,     &
92                                     free_cosp_sgradar      => destroy_cosp_sgradar,     &
93                                     free_cosp_radarstats   => destroy_cosp_radarstats,  &
94                                     free_cosp_sglidar      => destroy_cosp_sglidar,     &
95                                     free_cosp_lidarstats   => destroy_cosp_lidarstats,  &
96                                     free_cosp_isccp        => destroy_cosp_isccp,       &
97                                     free_cosp_misr         => destroy_cosp_misr,        &
98                                     free_cosp_rttov        => destroy_cosp_rttov,       &
99                                     free_cosp_modis        => destroy_cosp_modis,       &
100                                     free_cosp_vgrid        => destroy_cosp_vgrid,       &
101                                     free_cosp_subgrid      => destroy_cosp_subgrid,     &
102                                     construct_cosp_subgrid,cosp_config,cosp_subgrid,    &
103                                     cosp_sglidar,cosp_lidarstats,                       &
104                                     construct_cosp_lidarstats,construct_cosp_sglidar,   &
105                                     cosp_isccp,construct_cosp_isccp,cosp_misr,          &
106                                     construct_cosp_misr,cosp_rttov,construct_cosp_rttov,&
107                                     cosp_sgradar,cosp_radarstats,                       &
108                                     construct_cosp_radarstats,construct_cosp_sgradar,   &
109                                     cosp_modis,construct_cosp_modis,                    &
110                                     cosp_vgrid,I_CVCLIQ,I_LSCLIQ,I_CVCICE,I_LSCICE,     &
111                                     I_LSRAIN,I_LSSNOW,I_LSGRPL,I_CVRAIN,I_CVSNOW
112
113  use cosp_read_otputkeys
114
115  IMPLICIT NONE
116
117  ! Local variables
118  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
119  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
120  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
121  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
122  integer, save :: Npoints      ! Number of gridpoints
123!$OMP THREADPRIVATE(Npoints)
124  integer, save :: Nlevels      ! Number of levels
125  Integer :: Nptslmdz,Nlevlmdz ! Nb de points issus de physiq.F
126  integer, save :: Nlr          ! Number of levels in statistical outputs
127  integer, save :: Npoints_it   ! Max number of gridpoints to be processed in one iteration
128!  integer :: i
129  type(cosp_config),save :: cfg   ! Configuration options
130!$OMP THREADPRIVATE(cfg)
131  type(cosp_gridbox) :: gbx ! Gridbox information. Input for COSP
132  type(cosp_subgrid) :: sgx     ! Subgrid outputs
133  type(cosp_sgradar) :: sgradar ! Output from radar simulator
134  type(cosp_sglidar) :: sglidar ! Output from lidar simulator
135  type(cosp_isccp)   :: isccp   ! Output from ISCCP simulator
136!! AI rajout modis
137  type(cosp_modis)   :: modis   ! Output from MODIS simulator
138!!
139  type(cosp_misr)    :: misr    ! Output from MISR simulator
140!! AI rajout rttovs
141  type(cosp_rttov)   :: rttov   ! Output from RTTOV
142!!
143  type(cosp_vgrid)   :: vgrid   ! Information on vertical grid of stats
144  type(cosp_radarstats) :: stradar ! Summary statistics from radar simulator
145  type(cosp_lidarstats) :: stlidar ! Summary statistics from lidar simulator
146
147  integer :: t0,t1,count_rate,count_max
148  integer :: Nlon,Nlat
149  real(wp),save :: radar_freq,k2,ZenAng,co2,ch4,n2o,co,emsfc_lw
150 
151!$OMP THREADPRIVATE(emsfc_lw)
152  integer,dimension(RTTOV_MAX_CHANNELS),save :: Channels
153  real(wp),dimension(RTTOV_MAX_CHANNELS),save :: Surfem
154  integer, save :: surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay
155  integer, save :: Nprmts_max_hydro,Naero,Nprmts_max_aero,lidar_ice_type
156  integer, save :: platform,satellite,Instrument,Nchannels
157  logical, save :: use_vgrid,csat_vgrid,use_precipitation_fluxes,use_reff
158
159! Declaration necessaires pour les sorties IOIPSL
160  integer :: ii
161  real    :: ecrit_day,ecrit_hf,ecrit_mth, missing_val
162  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, ok_all_xml
163
164  logical, save :: debut_cosp=.true.
165!$OMP THREADPRIVATE(debut_cosp)
166
167  logical, save :: first_write=.true.
168!$OMP THREADPRIVATE(first_write)
169
170!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
171  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
172  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
173                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
174                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
175  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
176  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit         
177  real,dimension(Nlevlmdz)        :: presnivs
178  integer                         :: itap,k,ip
179  real                            :: dtime,freq_cosp
180  real,dimension(2)               :: time_bnds
181
182  double precision                            :: d_dtime
183  double precision,dimension(2)               :: d_time_bnds
184
185   namelist/COSP_INPUT/overlap,isccp_topheight,isccp_topheight_direction, &
186              npoints_it,ncolumns,use_vgrid,Nlr,csat_vgrid, &
187              radar_freq,surface_radar,use_mie_tables, &
188              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
189              lidar_ice_type,use_precipitation_fluxes,use_reff, &
190              platform,satellite,Instrument,Nchannels, &
191              Channels,Surfem,ZenAng,co2,ch4,n2o,co
192
193!---------------- End of declaration of variables --------------
194
195!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
196! Read namelist with COSP inputs
197!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
198
199 print*,'Entree phys_cosp2'
200 if (debut_cosp) then
201  NPoints=Nptslmdz
202  Nlevels=Nlevlmdz
203 
204! Lecture du namelist input
205!  CALL read_cosp_input
206   IF (is_master) THEN
207      OPEN(10,file=cosp_input_nl,status='old')
208      READ(10,nml=cosp_input)
209      CLOSE(10)
210   ENDIF
211!$OMP BARRIER
212    CALL bcast(overlap)
213    CALL bcast(isccp_topheight)
214    CALL bcast(isccp_topheight_direction)
215    CALL bcast(npoints_it)
216    CALL bcast(ncolumns)
217    CALL bcast(use_vgrid)
218    CALL bcast(Nlr)
219    CALL bcast(csat_vgrid)
220    CALL bcast(radar_freq)
221    CALL bcast(surface_radar)
222    CALL bcast(use_mie_tables)
223    CALL bcast(use_gas_abs)
224    CALL bcast(do_ray)
225    CALL bcast(melt_lay)
226    CALL bcast(k2)
227    CALL bcast(Nprmts_max_hydro)
228    CALL bcast(Naero)
229    CALL bcast(Nprmts_max_aero)
230    CALL bcast(lidar_ice_type)
231    CALL bcast(use_precipitation_fluxes)
232    CALL bcast(use_reff)
233    CALL bcast(platform)
234    CALL bcast(satellite)
235    CALL bcast(Instrument)
236    CALL bcast(Nchannels)
237    CALL bcast(Channels)
238    CALL bcast(Surfem)
239    CALL bcast(ZenAng)
240    CALL bcast(co2)
241    CALL bcast(ch4)
242    CALL bcast(n2o)
243    CALL bcast(co)
244
245    print*,'ok read  cosp_input_nl'
246
247! Clefs Outputs initialisation
248 IF (using_xios) THEN
249  call cosp_outputkeys_init(cfg)
250 ELSE
251   call read_cosp_output_nl(itap,cosp_output_nl,cfg)
252 ENDIF
253
254  print*,' Cles des differents simulateurs cosp a itap :',itap
255  print*,'Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lmodis_sim,Lrttov_sim,Lstats', &
256          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lmodis_sim, &
257          cfg%Lrttov_sim,cfg%Lstats
258
259    if (overlaplmdz.ne.overlap) then
260       print*,'Attention overlaplmdz different de overlap lu dans namelist '
261    endif
262   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
263
264  endif ! debut_cosp
265
266!!! Ici on modifie les cles logiques pour les outputs selon les champs actives dans les .xml
267  if ((itap.gt.1).and.(first_write))then
268   
269    IF (using_xios)   call read_xiosfieldactive(cfg)
270    first_write=.false.
271
272    print*,' Cles des differents simulateurs cosp a itap :',itap
273    print*,'Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lmodis_sim,Lrttov_sim,Lstats', &
274          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lmodis_sim, &
275          cfg%Lrttov_sim,cfg%Lstats
276  endif
277
278  time_bnds(1) = dtime-dtime/2.
279  time_bnds(2) = dtime+dtime/2.
280
281  d_time_bnds=time_bnds
282  d_dtime=dtime
283
284!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
285! Allocate memory for gridbox type
286!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
287! AI mars 2017
288!        print *, 'Allocating memory for gridbox type...'
289
290! Surafce emissivity
291        emsfc_lw = 1.
292
293        call construct_cosp_gridbox(d_dtime,d_time_bnds,radar_freq,surface_radar,use_mie_tables,use_gas_abs, &
294                                    do_ray,melt_lay,k2, &
295                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
296                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
297                                    use_precipitation_fluxes,use_reff, &
298                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
299                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
300       
301!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
302! Here code to populate input structure
303!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
304
305!        print *, 'Populating input structure...'
306        gbx%longitude = lon
307        gbx%latitude = lat
308
309        gbx%p = p !
310        gbx%ph = ph
311        gbx%zlev = phi/9.81
312
313        zlev_half(:,1) = phis(:)/9.81
314        do k = 2, Nlevels
315          do ip = 1, Npoints
316           zlev_half(ip,k) = phi(ip,k)/9.81 + &
317               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
318          enddo
319        enddo
320        gbx%zlev_half = zlev_half
321
322        gbx%T = T
323        gbx%q = rh*100.
324        gbx%sh = sh
325! On ne veut pas que cosp distingue les nuages stratiformes et convectifs
326! on passe les contenus totaux (conv+strat)
327        gbx%cca = 0. !convective_cloud_amount (1)
328        gbx%tca = tca ! total_cloud_amount (1)
329        gbx%psfc = ph(:,1) !pression de surface
330        gbx%skt  = skt !Skin temperature (K)
331
332        do ip = 1, Npoints
333          if (fracTerLic(ip).ge.0.5) then
334             gbx%land(ip) = 1.
335          else
336             gbx%land(ip) = 0.
337          endif
338        enddo
339        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
340! A voir l equivalent LMDZ (u10m et v10m)
341        gbx%u_wind  = u_wind !eastward_wind (m s-1)
342        gbx%v_wind  = v_wind !northward_wind
343
344! sunlit calcule a partir de la fraction d ensoleillement par jour
345!      do ip = 1, Npoints
346!        if (sunlit(ip).le.0.) then
347!           gbx%sunlit(ip)=0.
348!        else
349!           gbx%sunlit(ip)=1.
350!        endif
351!      enddo
352       gbx%sunlit=sunlit
353
354! A voir l equivalent LMDZ
355  mr_ccliq = 0.0
356  mr_ccice = 0.0
357        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
358        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
359        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
360        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
361! A revoir
362        fl_lsrain = fl_lsrainI + fl_ccrainI
363        fl_lssnow = fl_lssnowI + fl_ccsnowI
364        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
365        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
366!  A voir l equivalent LMDZ
367        fl_lsgrpl=0.
368        fl_ccsnow = 0.
369        fl_ccrain = 0.
370        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
371        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
372        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
373
374     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
375     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
376!! AI A revoir
377     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
378     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
379
380        ! ISCCP simulator
381        gbx%dtau_s   = dtau_s
382        gbx%dtau_c   = 0.
383        gbx%dem_s    = dem_s
384        gbx%dem_c   = 0.
385
386!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
387        ! Define new vertical grid
388!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
389!        print *, 'Defining new vertical grid...'
390        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
391        print*,'Nlvgrid,use_vgrid,csat_vgrid ',Nlr,use_vgrid,csat_vgrid
392        print*,'vgrid%z, vgrid%mz ',vgrid%z, vgrid%mz
393!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
394       ! Allocate memory for other types
395!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
396!        print *, 'Allocating memory for other types...'
397        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
398        if (cfg%Lradar_sim) call construct_cosp_sgradar(Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
399        if (cfg%Lradar_sim) call construct_cosp_radarstats(Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
400        if (cfg%Llidar_sim) call construct_cosp_sglidar(Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
401        if (cfg%Llidar_sim) call construct_cosp_lidarstats(Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
402        if (cfg%Lisccp_sim) call construct_cosp_isccp(Npoints,Ncolumns,Nlevels,isccp)
403!! AI rajout
404        if (cfg%Lmodis_sim) call construct_cosp_modis(Npoints,modis)
405!!
406        if (cfg%Lmisr_sim) call construct_cosp_misr(Npoints,misr)
407!        call construct_cosp_rttov(cfg,Npoints,Nchannels,rttov)
408
409    if (debut_cosp) then
410      !$OMP MASTER
411        print *, ' Open outpts files and define axis'
412        call cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp, &
413                              ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml, &
414                              ecrit_mth, ecrit_day, ecrit_hf, use_vgrid, vgrid)
415      !$OMP END MASTER
416      !$OMP BARRIER
417       debut_cosp=.false.
418      endif ! debut_cosp
419
420!!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
421        ! Call simulator
422!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
423      print*,'call simulateur'
424       call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,   &
425                                 isccp,misr,modis,rttov,stradar,stlidar)
426
427
428!!!!!!!!!!!!!!!!!! Ecreture des sorties Cosp !!!!!!!!!!!!!!r!!!!!!:!!!!!
429
430       print *, 'Calling write output'
431        call cosp_output_write(Nlevlmdz, Npoints, Ncolumns, itap, dtime, freq_COSP, missing_val, &
432                               cfg, gbx, vgrid, sglidar, sgradar, stlidar, stradar, &
433                               isccp, misr, modis)
434!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
435        ! Deallocate memory in derived types
436!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
437      print *, 'Deallocating memory...'
438     call free_cosp_gridbox(gbx)
439      print *, 'free_cosp_gridbox'
440     call free_cosp_subgrid(sgx)
441      print *, 'ffree_cosp_subgrid'
442     call free_cosp_vgrid(vgrid)
443      print *, 'free_cosp_vgrid'
444     if (cfg%Lradar_sim) call free_cosp_sgradar(sgradar)
445     if (cfg%Lradar_sim) call free_cosp_radarstats(stradar)
446     if (cfg%Llidar_sim) call free_cosp_sglidar(sglidar)
447     if (cfg%Llidar_sim) call free_cosp_lidarstats(stlidar)
448     if (cfg%Lisccp_sim) call free_cosp_isccp(isccp)
449     if (cfg%Lmisr_sim)  call free_cosp_misr(misr)
450     if (cfg%Lmodis_sim) call free_cosp_modis(modis)
451     if (cfg%Lrttov_sim) call free_cosp_rttov(rttov)
452     print *, 'End phys_cosp2'
453!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
454  ! Time in s. Only for testing purposes
455!  call system_clock(t1,count_rate,count_max)
456!  print *,(t1-t0)*1.0/count_rate
457 
458end subroutine phys_cosp2
Note: See TracBrowser for help on using the repository browser.