source: LMDZ4/branches/LMDZ4_AR5/libf/cosp/phys_cosp.F90 @ 1928

Last change on this file since 1928 was 1928, checked in by idelkadi, 10 years ago
  • Correction du bug dans l'interface avec la physique : On ne veut pas que la distinction entre les nuages convectifs et stratiformes soit 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
  • Rajouts d'autre variables diagnostiques : sunlit (=0 si nuit et =1 si jour) parasol_crefl (reflectance integree)
File size: 20.0 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 
79  IMPLICIT NONE
80
81  ! Local variables
82  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
83  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
84  character(len=512), save :: finput ! Input file name
85  character(len=512), save :: cmor_nl
86  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
87  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
88!  integer,parameter :: Ncollmdz=20
89  integer,parameter :: Ncolmax=100
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,idayref
122  real    :: zjulian,zstoday,zstomth,zstohf,zout,ecrit_day,ecrit_hf,ecrit_mth
123  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP
124  integer :: nhori,nvert,nvertp,nvertisccp,nvertm,nvertcol
125  integer, save :: nid_day_cosp,nid_mth_cosp,nid_hf_cosp
126!$OMP THREADPRIVATE(nid_day_cosp,nid_mth_cosp,nid_hf_cosp)
127  logical, save :: debut_cosp=.true.
128!$OMP THREADPRIVATE(debut_cosp)
129  integer :: itau_wcosp
130  character(len=2) :: str2
131  real,dimension(Ncolmax) :: column_ax
132  character(len=10),save,dimension(Ncolmax) :: chcol
133
134  integer, save :: Nlevout
135!$OMP THREADPRIVATE(Nlevout)
136
137  include "dimensions.h"
138  include "temps.h" 
139 
140!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
141  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
142  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
143                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
144                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
145  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
146  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit   
147  real,dimension(Nlevlmdz)        :: presnivs
148  integer                         :: itap,k,ip
149  real                            :: dtime,freq_cosp
150  logical, parameter              :: lCOSP=.FALSE.
151
152  real, dimension(Nptslmdz,PARASOL_NREFL) :: parasolcrefl, Ncref
153
154!
155   namelist/COSP_INPUT/cmor_nl,overlap,isccp_topheight,isccp_topheight_direction, &
156              npoints_it,ncolumns,nlevels,use_vgrid,nlr,csat_vgrid,finput, &
157              radar_freq,surface_radar,use_mie_tables, &
158              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
159              lidar_ice_type,use_precipitation_fluxes,use_reff, &
160              platform,satellite,Instrument,Nchannels, &
161              Channels,Surfem,ZenAng,co2,ch4,n2o,co
162
163!---------------- End of declaration of variables --------------
164   
165
166!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167! Read namelist with COSP inputs
168!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169
170 if (debut_cosp) then
171  NPoints=Nptslmdz
172  Nlevels=Nlevlmdz
173 
174! Lecture du namelist input
175  CALL read_cosp_input
176
177  do ii=1,Ncolumns
178    write(str2,'(i2.2)')ii
179    chcol(ii)="c"//str2
180    column_ax(ii) = real(ii)
181  enddo
182
183! Clefs Outputs
184  call read_cosp_output_nl(cosp_output_nl,cfg)
185
186    if (overlaplmdz.ne.overlap) then
187       print*,'Attention overlaplmdz different de overlap lu dans namelist '
188    endif
189   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
190
191  print*,' Cles sorties cosp :'
192  print*,' Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lrttov_sim', &
193          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lrttov_sim
194
195  endif ! debut_cosp
196
197  print*,'Debut phys_cosp itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf ', &
198          itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf
199!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
200! Allocate local arrays
201!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
202!        call system_clock(t0,count_rate,count_max) !!! Only for testing purposes
203       
204       
205!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
206! Allocate memory for gridbox type
207!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
208        print *, 'Allocating memory for gridbox type...'
209
210        call construct_cosp_gridbox(dble(itap),radar_freq,surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay,k2, &
211                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
212                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
213                                    use_precipitation_fluxes,use_reff, &
214                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
215                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
216       
217!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
218! Here code to populate input structure
219!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
220
221        print *, 'Populating input structure...'
222        gbx%longitude = lon
223        gbx%latitude = lat
224
225        gbx%p = p !
226        gbx%ph = ph
227        gbx%zlev = phi/9.81
228
229        zlev_half(:,1) = phis(:)/9.81
230        do k = 2, Nlevels
231          do ip = 1, Npoints
232           zlev_half(ip,k) = phi(ip,k)/9.81 + &
233               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
234          enddo
235        enddo
236        gbx%zlev_half = zlev_half
237
238        gbx%T = T
239        gbx%q = rh*100.
240        gbx%sh = sh
241        gbx%cca = 0.
242        gbx%tca = tca ! total_cloud_amount (1)
243        gbx%psfc = ph(:,1) !pression de surface
244        gbx%skt  = skt !Skin temperature (K)
245
246        do ip = 1, Npoints
247          if (fracTerLic(ip).ge.0.5) then
248             gbx%land(ip) = 1.
249          else
250             gbx%land(ip) = 0.
251          endif
252        enddo
253        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
254! A voir l equivalent LMDZ (u10m et v10m)
255        gbx%u_wind  = u_wind !eastward_wind (m s-1)
256        gbx%v_wind  = v_wind !northward_wind
257
258      do ip = 1, Npoints
259        if (sunlit(ip).le.0.) then
260           gbx%sunlit(ip)=0
261        else
262           gbx%sunlit(ip)=1
263        endif
264      enddo
265
266! A voir l equivalent LMDZ
267  mr_ccliq = 0.0
268  mr_ccice = 0.0
269        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
270        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
271        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
272        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
273! A revoir
274        fl_lsrain = fl_lsrainI + fl_ccrainI
275        fl_lssnow = fl_lssnowI + fl_ccsnowI
276        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
277        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
278!  A voir l equivalent LMDZ
279        fl_lsgrpl=0.
280        fl_ccsnow = 0.
281        fl_ccrain = 0.
282        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
283        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
284        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
285
286     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
287     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
288     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
289     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
290        ! ISCCP simulator
291        gbx%dtau_s   = dtau_s
292        gbx%dtau_c   = 0.
293        gbx%dem_s    = dem_s
294        gbx%dem_c    = 0.
295
296! Surafce emissivity
297       emsfc_lw = 1.
298               
299!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
300        ! Define new vertical grid
301!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
302        print *, 'Defining new vertical grid...'
303        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
304
305 if (debut_cosp) then
306! Creer le fichier de sorie, definir les variable de sortie
307  ! Axe verticale (Pa ou Km)
308     Nlevout = vgrid%Nlvgrid
309   
310        do ii=1,Ncolumns
311          column_ax(ii) = real(ii)
312        enddo
313
314 if (ok_mensuelCOSP) then
315     include "ini_histmthCOSP.h"
316 endif
317 if (ok_journeCOSP) then
318     include "ini_histdayCOSP.h"
319 endif
320 if (ok_hfCOSP) then
321     include "ini_histhfCOSP.h"
322 endif
323
324   debut_cosp=.false.
325  endif ! debut_cosp
326
327!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
328       ! Allocate memory for other types
329!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
330        print *, 'Allocating memory for other types...'
331        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
332        call construct_cosp_sgradar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
333        call construct_cosp_radarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
334        call construct_cosp_sglidar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
335        call construct_cosp_lidarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
336        call construct_cosp_isccp(cfg,Npoints,Ncolumns,Nlevels,isccp)
337        call construct_cosp_misr(cfg,Npoints,misr)
338!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
339        ! Call simulator
340!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
341        print *, 'Calling simulator...'
342        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,stradar,stlidar)
343!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
344        ! Write outputs to CMOR-compliant NetCDF
345!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
346
347! A traiter le cas ou l on a des valeurs indefinies
348! Attention teste
349
350! if(1.eq.0)then
351
352
353   do k = 1,Nlevout
354     do ip = 1,Npoints
355     if(stlidar%lidarcld(ip,k).eq.R_UNDEF)then
356      stlidar%lidarcld(ip,k)=0.
357     endif
358     enddo
359
360
361     do ii= 1,SR_BINS
362      do ip = 1,Npoints
363       if(stlidar%cfad_sr(ip,ii,k).eq.R_UNDEF)then
364        stlidar%cfad_sr(ip,ii,k)=0.
365       endif
366      enddo
367     enddo
368   enddo   
369   
370  do ip = 1,Npoints
371   do k = 1,Nlevlmdz
372     if(sglidar%beta_mol(ip,k).eq.R_UNDEF)then
373      sglidar%beta_mol(ip,k)=0.
374     endif
375   
376     do ii= 1,Ncolumns
377       if(sglidar%beta_tot(ip,ii,k).eq.R_UNDEF)then
378        sglidar%beta_tot(ip,ii,k)=0.
379       endif 
380     enddo
381
382    enddo    !k = 1,Nlevlmdz
383   enddo     !ip = 1,Npoints
384
385   do k = 1,LIDAR_NCAT
386    do ip = 1,Npoints
387     if(stlidar%cldlayer(ip,k).eq.R_UNDEF)then
388      stlidar%cldlayer(ip,k)=0.
389     endif
390    enddo
391   enddo
392
393! endif
394
395   do ip = 1,Npoints
396    if(isccp%totalcldarea(ip).eq.-1.E+30)then
397      isccp%totalcldarea(ip)=0.
398    endif
399    if(isccp%meanptop(ip).eq.-1.E+30)then
400      isccp%meanptop(ip)=0.
401    endif
402    if(isccp%meantaucld(ip).eq.-1.E+30)then
403      isccp%meantaucld(ip)=0.
404    endif
405    if(isccp%meanalbedocld(ip).eq.-1.E+30)then
406      isccp%meanalbedocld(ip)=0.
407    endif
408    if(isccp%meantb(ip).eq.-1.E+30)then
409      isccp%meantb(ip)=0.
410    endif
411    if(isccp%meantbclr(ip).eq.-1.E+30)then
412      isccp%meantbclr(ip)=0.
413    endif
414
415    do k=1,7
416     do ii=1,7
417     if(isccp%fq_isccp(ip,ii,k).eq.-1.E+30)then
418      isccp%fq_isccp(ip,ii,k)=0.
419     endif
420     enddo
421    enddo
422
423    do ii=1,Ncolumns
424     if(isccp%boxtau(ip,ii).eq.-1.E+30)then
425       isccp%boxtau(ip,ii)=0.
426     endif
427    enddo
428
429    do ii=1,Ncolumns
430     if(isccp%boxptop(ip,ii).eq.-1.E+30)then
431       isccp%boxptop(ip,ii)=0.
432     endif
433    enddo
434   enddo
435
436   do k=1,PARASOL_NREFL
437    do ii=1, Npoints
438     if (stlidar%cldlayer(ii,4).gt.0.01) then
439        parasolcrefl(ii,k)=(stlidar%parasolrefl(ii,k)-0.03*(1.-stlidar%cldlayer(ii,4)))/ &
440                   stlidar%cldlayer(ii,4)
441        Ncref(ii,k) = 1.
442     else
443        parasolcrefl(ii,k)=0.
444        Ncref(ii,k) = 0.
445     endif
446    enddo
447   enddo
448
449 if (ok_mensuelCOSP) then
450  include "write_histmthCOSP.h"
451 endif
452 if (ok_journeCOSP) then
453  include "write_histdayCOSP.h"
454 endif
455 if (ok_hfCOSP ) then
456  include "write_histhfCOSP.h"
457 endif
458
459!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
460        ! Deallocate memory in derived types
461!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
462        print *, 'Deallocating memory...'
463        call free_cosp_gridbox(gbx)
464        call free_cosp_subgrid(sgx)
465        call free_cosp_sgradar(sgradar)
466        call free_cosp_radarstats(stradar)
467        call free_cosp_sglidar(sglidar)
468        call free_cosp_lidarstats(stlidar)
469        call free_cosp_isccp(isccp)
470        call free_cosp_misr(misr)
471        call free_cosp_vgrid(vgrid) 
472 
473!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
474  ! Time in s. Only for testing purposes
475!  call system_clock(t1,count_rate,count_max)
476!  print *,(t1-t0)*1.0/count_rate
477 
478  CONTAINS
479 
480  SUBROUTINE read_cosp_input
481   
482    IF (is_master) THEN
483      OPEN(10,file=cosp_input_nl,status='old')
484      READ(10,nml=cosp_input)
485      CLOSE(10)
486    ENDIF
487    CALL bcast(cmor_nl)
488    CALL bcast(overlap)
489    CALL bcast(isccp_topheight)
490    CALL bcast(isccp_topheight_direction)
491    CALL bcast(npoints_it)
492    CALL bcast(ncolumns)
493    CALL bcast(nlevels)
494    CALL bcast(use_vgrid)
495    CALL bcast(nlr)
496    CALL bcast(csat_vgrid)
497    CALL bcast(finput)
498    CALL bcast(radar_freq)
499    CALL bcast(surface_radar)
500    CALL bcast(use_mie_tables)
501    CALL bcast(use_gas_abs)
502    CALL bcast(do_ray)
503    CALL bcast(melt_lay)
504    CALL bcast(k2)
505    CALL bcast(Nprmts_max_hydro)
506    CALL bcast(Naero)
507    CALL bcast(Nprmts_max_aero)
508    CALL bcast(lidar_ice_type)
509    CALL bcast(use_precipitation_fluxes)
510    CALL bcast(use_reff)
511    CALL bcast(platform)
512    CALL bcast(satellite)
513    CALL bcast(Instrument)
514    CALL bcast(Nchannels)
515    CALL bcast(Channels)
516    CALL bcast(Surfem)
517    CALL bcast(ZenAng)
518    CALL bcast(co2)
519    CALL bcast(ch4)
520    CALL bcast(n2o)
521    CALL bcast(co)
522!$OMP BARRIER 
523  END SUBROUTINE read_cosp_input
524
525end subroutine phys_cosp
Note: See TracBrowser for help on using the repository browser.