source: LMDZ6/branches/IPSLCM6.0.15/libf/phylmd/cosp/phys_cosp.F90 @ 3246

Last change on this file since 3246 was 3246, checked in by idelkadi, 6 years ago

Re-ecriture de la routine de gestion des cles de sorties Cosp :

ecriture d'un module cosp_read_otputkeys.F90 incluant 3 routines :

  • une subroutine pour initialiser au 1er appel a Cosp les Cles
  • une subroutine appelee au 2e passage dans Cosp, permettant de piloter les cles de sorties en fonction de ce qui est demande dans les fichiers xml (dans le cas ou les sorties sont geres par Xios)
  • une routine appelee au 2e passage dans Cosp, permettant de lire les cles dans le ficher cosp_output_nl.txt dans le cas ou les sorties ne sont pas geres par Xios

Nettoyage dans le module cosp_output_write.F90

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