source: LMDZ4/trunk/libf/cosp/phys_cosp.F90 @ 1368

Last change on this file since 1368 was 1368, checked in by idelkadi, 14 years ago
  • Rajout des cles logiques pour activer les sorties COSP par frequence
  • Reglages des niveaux de sorties pour les variables CFMIP
File size: 20.8 KB
Line 
1! Simulateur COSP : Cfmip Observation Simulator Package
2! ISCCP, Radar (QuickBeam), Lidar et Parasol (ACTSIM), MISR, RTTOVS
3!Idelkadi Abderrahmane Aout-Septembre 2009
4
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, &
10                        ref_liq,ref_ice,fracTerLic,u_wind,v_wind,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! ph,                                   !pression pour chaque inter-couche
28! p,                                    !Pression aux milieux des couches
29! skt,t,                                !Temp au sol et temp 3D
30! sh,                                   !Humidite specifique
31! rh,                                   !Humidite relatif
32! tca,                                  !Fraction nuageuse
33! cca                                   !Fraction nuageuse convective
34! mr_lsliq,                             !Liq Cloud water content
35! mr_lsice,                             !Ice Cloud water content
36! mr_ccliq,                             !Convective Cloud Liquid water content 
37! mr_ccice,                             !Cloud ice water content
38! fl_lsrain,                            !Large scale precipitation lic
39! fl_lssnow,                            !Large scale precipitation ice
40! fl_ccrain,                            !Convective precipitation lic
41! fl_ccsnow,                            !Convective precipitation ice
42! mr_ozone,                             !Concentration ozone (Kg/Kg)
43! dem_s                                 !Cloud optical emissivity
44! dtau_s                                !Cloud optical thickness
45! emsfc_lw = 1.                         !Surface emissivity dans radlwsw.F90
46
47!!! Outputs :
48! calipso2D,                            !Lidar Low/heigh/Mean/Total-level Cloud Fraction
49! calipso3D,                            !Lidar Cloud Fraction (532 nm)
50! cfadlidar,                            !Lidar Scattering Ratio CFAD (532 nm)
51! parasolrefl,                          !PARASOL-like mono-directional reflectance
52! atb,                                  !Lidar Attenuated Total Backscatter (532 nm)
53! betamol,                              !Lidar Molecular Backscatter (532 nm)
54! cfaddbze,                             !Radar Reflectivity Factor CFAD (94 GHz)
55! clcalipso2,                           !Cloud frequency of occurrence as seen by CALIPSO but not CloudSat
56! dbze,                                 !Efective_reflectivity_factor
57! cltlidarradar,                        !Lidar and Radar Total Cloud Fraction
58! clMISR,                               !Cloud Fraction as Calculated by the MISR Simulator
59! clisccp2,                             !Cloud Fraction as Calculated by the ISCCP Simulator
60! boxtauisccp,                          !Optical Depth in Each Column as Calculated by the ISCCP Simulator
61! boxptopisccp,                         !Cloud Top Pressure in Each Column as Calculated by the ISCCP Simulator
62! tclisccp,                             !Total Cloud Fraction as Calculated by the ISCCP Simulator
63! ctpisccp,                             !Mean Cloud Top Pressure as Calculated by the ISCCP Simulator
64! tauisccp,                             !Mean Optical Depth as Calculated by the ISCCP Simulator
65! albisccp,                             !Mean Cloud Albedo as Calculated by the ISCCP Simulator
66! meantbisccp,                          !Mean all-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
67! meantbclrisccp                        !Mean clear-sky 10.5 micron brightness temperature as calculated by the ISCCP Simulator
68!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
69
70  USE MOD_COSP_CONSTANTS
71  USE MOD_COSP_TYPES
72  USE MOD_COSP
73  USE mod_phys_lmdz_para
74  USE mod_grid_phy_lmdz
75  use ioipsl
76  use iophy
77 
78  IMPLICIT NONE
79
80  ! Local variables
81  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
82  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
83  character(len=512), save :: finput ! Input file name
84  character(len=512), save :: cmor_nl
85  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
86  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
87  integer,parameter :: Ncollmdz=20
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,idayref
120  real    :: zjulian,zstoday,zstomth,zstohf,zout,ecrit_day,ecrit_hf,ecrit_mth
121  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP
122  integer :: nhori,nvert,nvertp,nvertisccp,nvertm,nvertcol
123  integer, save :: nid_day_cosp,nid_mth_cosp,nid_hf_cosp
124!$OMP THREADPRIVATE(nid_day_cosp,nid_mth_cosp,nid_hf_cosp)
125  logical, save :: debut_cosp=.true.
126!$OMP THREADPRIVATE(debut_cosp)
127  integer :: itau_wcosp
128  character(len=10),dimension(Ncollmdz),parameter :: chcol=(/'c01','c02','c03','c04','c05','c06','c07','c08','c09','c10', &
129                                                   'c11','c12','c13','c14','c15','c16','c17','c18','c19','c20'/)
130  real,dimension(Ncollmdz) :: column_ax
131  integer, save :: Nlevout
132!$OMP THREADPRIVATE(Nlevout)
133
134  include "dimensions.h"
135  include "temps.h" 
136 
137!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
138  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
139!  real,dimension(Npoints,Nlevels) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
140  real,dimension(Nptslmdz,Nlevlmdz) :: height,phi,p,ph,T,sh,rh,tca,cca,mr_lsliq,mr_lsice,mr_ccliq,mr_ccice, &
141                                     fl_lsrain,fl_lssnow,fl_ccrain,fl_ccsnow,fl_lsgrpl, &
142                                     zlev,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
143  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
144  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind
145  real,dimension(Nlevlmdz)        :: presnivs
146  integer                         :: itap,k,ip
147  real                            :: dtime,freq_cosp
148
149!
150   namelist/COSP_INPUT/cmor_nl,overlap,isccp_topheight,isccp_topheight_direction, &
151              npoints_it,ncolumns,nlevels,use_vgrid,nlr,csat_vgrid,finput, &
152              radar_freq,surface_radar,use_mie_tables, &
153              use_gas_abs,do_ray,melt_lay,k2,Nprmts_max_hydro,Naero,Nprmts_max_aero, &
154              lidar_ice_type,use_precipitation_fluxes,use_reff, &
155              platform,satellite,Instrument,Nchannels, &
156              Channels,Surfem,ZenAng,co2,ch4,n2o,co
157
158!---------------- End of declaration of variables --------------
159   
160
161!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162! Read namelist with COSP inputs
163!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164
165 if (debut_cosp) then
166  NPoints=Nptslmdz
167! Lecture du namelist input
168  CALL read_cosp_input
169
170! Clefs Outputs
171  call read_cosp_output_nl(cosp_output_nl,cfg)
172
173    if ( (Ncollmdz.ne.Ncolumns).or. (Nlevlmdz.ne.Nlevels) ) then
174       print*,'Nb points Horiz, Vert, Sub-col passes par physiq.F = ', &
175               Nptslmdz, Nlevlmdz, Ncollmdz
176       print*,'Nb points Horiz, Vert, Sub-col lus dans namelist = ', &
177               Npoints, Nlevels, Ncolumns
178       print*,'Nb points Horiz, Vert, Sub-col passes par physiq.F est different de celui lu par namelist '
179       call abort
180    endif
181   
182    if (overlaplmdz.ne.overlap) then
183       print*,'Attention overlaplmdz different de overlap lu dans namelist '
184    endif
185   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
186
187  print*,' Cles sorties cosp :'
188  print*,' Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lrttov_sim', &
189          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lrttov_sim
190
191  endif ! debut_cosp
192
193  print*,'Debut phys_cosp itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf ', &
194          itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf
195!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
196! Allocate local arrays
197!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
198!        call system_clock(t0,count_rate,count_max) !!! Only for testing purposes
199       
200       
201!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
202! Allocate memory for gridbox type
203!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
204        print *, 'Allocating memory for gridbox type...'
205
206        call construct_cosp_gridbox(dble(itap),radar_freq,surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay,k2, &
207                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
208                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
209                                    use_precipitation_fluxes,use_reff, &
210                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
211                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
212       
213!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
214! Here code to populate input structure
215!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
216
217        print *, 'Populating input structure...'
218        gbx%longitude = lon
219        gbx%latitude = lat
220
221        gbx%p = p !
222        gbx%ph = ph
223        gbx%zlev_half = phi/9.81
224
225       do k = 1, Nlevels-1
226       do ip = 1, Npoints
227        zlev(ip,k) = phi(ip,k)/9.81 + (phi(ip,k+1)-phi(ip,k))/9.81 * (ph(ip,k)-ph(ip,k+1))/p(ip,k)
228       enddo
229       enddo
230       do ip = 1, Npoints
231        zlev(ip,Nlevels) = zlev(ip,Nlevels-1)+ 2.*(phi(ip,Nlevels)/9.81-zlev(ip,Nlevels-1))
232       END DO
233        gbx%zlev = zlev
234
235        gbx%T = T
236        gbx%q = rh*100.
237        gbx%sh = sh
238        gbx%cca = cca !convective_cloud_amount (1)
239        gbx%tca = tca ! total_cloud_amount (1)
240        gbx%psfc = ph(:,1) !pression de surface
241        gbx%skt  = skt !Skin temperature (K)
242
243        do ip = 1, Npoints
244          if (fracTerLic(ip).ge.0.5) then
245             gbx%land(ip) = 1.
246          else
247             gbx%land(ip) = 0.
248          endif
249        enddo
250        gbx%mr_ozone  = mr_ozone !mass_fraction_of_ozone_in_air (kg/kg)
251! A voir l equivalent LMDZ (u10m et v10m)
252        gbx%u_wind  = u_wind !eastward_wind (m s-1)
253        gbx%v_wind  = v_wind !northward_wind
254! Attention
255        gbx%sunlit  = 1
256
257! A voir l equivalent LMDZ
258  mr_ccliq = 0.0
259  mr_ccice = 0.0
260        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
261        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
262        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
263        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
264! A revoir
265        fl_lsrain = fl_lsrainI + fl_ccrainI
266        fl_lssnow = fl_lssnowI + fl_ccsnowI
267        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
268        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
269!  A voir l equivalent LMDZ
270        fl_lsgrpl=0.
271        fl_ccsnow = 0.
272        fl_ccrain = 0.
273        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
274        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
275        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
276
277!Attention Teste
278!       do k = 1, Nlevels
279!        do ip = 1, Npoints
280!!     liquid particles :
281!         radliq(ip,k) = 12.0e-06
282!         if (k.le.3) radliq(ip,k) = 11.0e-06
283
284!    ice particles :
285!        if ( (t(ip,k)-273.15).gt.-81.4 ) then
286!          radice(ip,k) = (0.71*(t(ip,k)-273.15)+61.29)*1e-6
287!        else
288!          radice(ip,k) = 3.5*1e-6
289!        endif
290!       END DO
291!      END DO
292
293!      gbx%Reff(:,:,I_LSCLIQ) = radliq
294!      gbx%Reff(:,:,I_LSCICE) = radice
295!      gbx%Reff(:,:,I_CVCLIQ) = radliq
296!      gbx%Reff(:,:,I_CVCICE) = radice
297!      print*,'radliq(1,:)=',radliq(1,:)
298!      print*,'radice(1,:)=',radice(1,:)
299
300     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
301     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
302     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
303     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
304!     print*,'ref_liq(1,:)=',ref_liq(1,:)*1e-6
305!     print*,'ref_liq(1,:)=',ref_ice(1,:)*1e-6
306
307        ! ISCCP simulator
308        gbx%dtau_s   = dtau_s
309        gbx%dtau_c   = 0.
310        gbx%dem_s    = dem_s
311        gbx%dem_c    = 0.
312
313! Surafce emissivity
314       emsfc_lw = 1.
315               
316!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
317        ! Define new vertical grid
318!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
319        print *, 'Defining new vertical grid...'
320        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
321
322 if (debut_cosp) then
323! Creer le fichier de sorie, definir les variable de sortie
324  ! Axe verticale (Pa ou Km)
325     Nlevout = vgrid%Nlvgrid
326   
327        do ii=1,Ncolumns
328          column_ax(ii) = real(ii)
329        enddo
330
331 if (ok_mensuelCOSP) then
332     include "ini_histmthCOSP.h"
333 endif
334 if (ok_journeCOSP) then
335     include "ini_histdayCOSP.h"
336 endif
337 if (ok_hfCOSP) then
338     include "ini_histhfCOSP.h"
339 endif
340
341!   print*,'Fin Initialisation des sorties COSP, debut_cosp =',debut_cosp
342!   print*,'R_UNDEF=',R_UNDEF
343
344   debut_cosp=.false.
345  endif ! debut_cosp
346
347!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
348       ! Allocate memory for other types
349!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
350        print *, 'Allocating memory for other types...'
351        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
352        call construct_cosp_sgradar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
353        call construct_cosp_radarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
354        call construct_cosp_sglidar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
355        call construct_cosp_lidarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
356        call construct_cosp_isccp(cfg,Npoints,Ncolumns,Nlevels,isccp)
357        call construct_cosp_misr(cfg,Npoints,misr)
358!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
359        ! Call simulator
360!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
361        print *, 'Calling simulator...'
362        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,stradar,stlidar)
363!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
364        ! Write outputs to CMOR-compliant NetCDF
365!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
366
367! A traiter le cas ou l on a des valeurs indefinies
368! Attention teste
369
370! if(1.eq.0)then
371
372
373   do k = 1,Nlevout
374     do ip = 1,Npoints
375     if(stlidar%lidarcld(ip,k).eq.R_UNDEF)then
376      stlidar%lidarcld(ip,k)=0.
377     endif
378     enddo
379
380
381     do ii= 1,SR_BINS
382      do ip = 1,Npoints
383       if(stlidar%cfad_sr(ip,ii,k).eq.R_UNDEF)then
384        stlidar%cfad_sr(ip,ii,k)=0.
385       endif
386      enddo
387     enddo
388   enddo   
389   
390  do ip = 1,Npoints
391   do k = 1,Nlevlmdz
392     if(sglidar%beta_mol(ip,k).eq.R_UNDEF)then
393      sglidar%beta_mol(ip,k)=0.
394     endif
395   
396     do ii= 1,Ncolumns
397       if(sglidar%beta_tot(ip,ii,k).eq.R_UNDEF)then
398        sglidar%beta_tot(ip,ii,k)=0.
399       endif 
400     enddo
401
402    enddo    !k = 1,Nlevlmdz
403   enddo     !ip = 1,Npoints
404
405   do k = 1,LIDAR_NCAT
406    do ip = 1,Npoints
407     if(stlidar%cldlayer(ip,k).eq.R_UNDEF)then
408      stlidar%cldlayer(ip,k)=0.
409     endif
410    enddo
411   enddo
412
413! endif
414
415   do ip = 1,Npoints
416    if(isccp%totalcldarea(ip).eq.-1.E+30)then
417      isccp%totalcldarea(ip)=0.
418    endif
419    if(isccp%meanptop(ip).eq.-1.E+30)then
420      isccp%meanptop(ip)=0.
421    endif
422    if(isccp%meantaucld(ip).eq.-1.E+30)then
423      isccp%meantaucld(ip)=0.
424    endif
425    if(isccp%meanalbedocld(ip).eq.-1.E+30)then
426      isccp%meanalbedocld(ip)=0.
427    endif
428    if(isccp%meantb(ip).eq.-1.E+30)then
429      isccp%meantb(ip)=0.
430    endif
431    if(isccp%meantbclr(ip).eq.-1.E+30)then
432      isccp%meantbclr(ip)=0.
433    endif
434
435    do k=1,7
436     do ii=1,7
437     if(isccp%fq_isccp(ip,ii,k).eq.-1.E+30)then
438      isccp%fq_isccp(ip,ii,k)=0.
439     endif
440     enddo
441    enddo
442
443    do ii=1,Ncolumns
444     if(isccp%boxtau(ip,ii).eq.-1.E+30)then
445       isccp%boxtau(ip,ii)=0.
446     endif
447    enddo
448
449    do ii=1,Ncolumns
450     if(isccp%boxptop(ip,ii).eq.-1.E+30)then
451       isccp%boxptop(ip,ii)=0.
452     endif
453    enddo
454   enddo
455
456 if (ok_mensuelCOSP) then
457  include "write_histmthCOSP.h"
458 endif
459 if (ok_journeCOSP) then
460  include "write_histdayCOSP.h"
461 endif
462 if (ok_hfCOSP ) then
463  include "write_histhfCOSP.h"
464 endif
465
466!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
467        ! Deallocate memory in derived types
468!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
469        print *, 'Deallocating memory...'
470        call free_cosp_gridbox(gbx)
471        call free_cosp_subgrid(sgx)
472        call free_cosp_sgradar(sgradar)
473        call free_cosp_radarstats(stradar)
474        call free_cosp_sglidar(sglidar)
475        call free_cosp_lidarstats(stlidar)
476        call free_cosp_isccp(isccp)
477        call free_cosp_misr(misr)
478        call free_cosp_vgrid(vgrid) 
479 
480!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
481  ! Time in s. Only for testing purposes
482!  call system_clock(t1,count_rate,count_max)
483!  print *,(t1-t0)*1.0/count_rate
484 
485  CONTAINS
486 
487  SUBROUTINE read_cosp_input
488   
489    IF (is_master) THEN
490      OPEN(10,file=cosp_input_nl,status='old')
491      READ(10,nml=cosp_input)
492      CLOSE(10)
493    ENDIF
494    CALL bcast(cmor_nl)
495    CALL bcast(overlap)
496    CALL bcast(isccp_topheight)
497    CALL bcast(isccp_topheight_direction)
498    CALL bcast(npoints_it)
499    CALL bcast(ncolumns)
500    CALL bcast(nlevels)
501    CALL bcast(use_vgrid)
502    CALL bcast(nlr)
503    CALL bcast(csat_vgrid)
504    CALL bcast(finput)
505    CALL bcast(radar_freq)
506    CALL bcast(surface_radar)
507    CALL bcast(use_mie_tables)
508    CALL bcast(use_gas_abs)
509    CALL bcast(do_ray)
510    CALL bcast(melt_lay)
511    CALL bcast(k2)
512    CALL bcast(Nprmts_max_hydro)
513    CALL bcast(Naero)
514    CALL bcast(Nprmts_max_aero)
515    CALL bcast(lidar_ice_type)
516    CALL bcast(use_precipitation_fluxes)
517    CALL bcast(use_reff)
518    CALL bcast(platform)
519    CALL bcast(satellite)
520    CALL bcast(Instrument)
521    CALL bcast(Nchannels)
522    CALL bcast(Channels)
523    CALL bcast(Surfem)
524    CALL bcast(ZenAng)
525    CALL bcast(co2)
526    CALL bcast(ch4)
527    CALL bcast(n2o)
528    CALL bcast(co)
529!$OMP BARRIER 
530  END SUBROUTINE read_cosp_input
531
532end subroutine phys_cosp
Note: See TracBrowser for help on using the repository browser.