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

Last change on this file since 1535 was 1535, checked in by musat, 13 years ago

Ajout flag LOGICAL lCOSP necessaire pour sortir un fichier stations
IM

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