source: LMDZ4/branches/LMDZ4V5.0-dev/libf/cosp/phys_cosp.F90 @ 1293

Last change on this file since 1293 was 1293, checked in by lguez, 14 years ago

1) Replaced calls to "float" by calls to "real" and "dble" in
"phys_cosp". "call construct_cosp_gridbox(float(itap),..." was a bug
since the corresponding dummy argument has the type "double
precision". (And "float" is not in the Fortran standard.)

2) Modifications for the program "create_etat0_limit" only, "gcm" is
not touched:

2.1) Removed generic interface "startget" in module "startvar". Replaced
calls to "startget" by calls to "startget_phys2d", "startget_dyn" and
"startget_phys1d".

2.2) Simplified the interface of "startget_dyn" and made it more secure by
removing arguments for sizes of arrays and using assumed-shape array
arguments.

2.3) Corrected bug in call to "startget_dyn" for northward velocity. See
ticket #25 in Trac.

2.4) Collected procedures "inter_barxy", "inter_barx", "inter_bary",
"ord_coord" and "ord_coordm" in a module. "inter_barxy" is the only
public procedure in the module. Rewrote those five procedures: removed
arguments for sizes of arrays; used assumed-shape array arguments;
used array expressions; translated some spaghetti code in "inter_bary"
into structured code; corrected bug in "inter_bary", described by
ticket #26 in Trac.

2.5) Corrected a bug in "inter_bary": "y0" was initialized at 0
instead of -90. This bug made values at the south pole wrong. (This bug
has been acknowledeged by P. Le Van.)

2.6) Made the variable "champ" in procedure "limit_netcdf" an array of
rank 2. Thus it can be an argument in the call to the newly-secured
"inter_barxy".

2.7) The files "start.nc", "startphy.nc" and "limit.nc" are impacted
by this revision. There is a significant change in the variable "vcov"
of "start.nc". See
http://web.lmd.jussieu.fr/~lglmd/Rev_1293/index.html. Note also the
changed value of "vcov" at the south pole. There is very little change
in other variables, in all three files.

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