source: LMDZ5/trunk/libf/phylmd/cosp/phys_cosp.F90

Last change on this file was 2835, checked in by idelkadi, 7 years ago

Supression des messages de sorties (print*) dans Cosp!

  • 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: 18.8 KB
RevLine 
[1262]1! Simulateur COSP : Cfmip Observation Simulator Package
[1526]2
[1262]3! ISCCP, Radar (QuickBeam), Lidar et Parasol (ACTSIM), MISR, RTTOVS
[2428]4!Idelkadi Abderrahmane Aout-Septembre 2009 First Version
5!Idelkadi Abderrahmane Nov 2015 version v1.4.0
[1262]6
[1368]7  subroutine phys_cosp( itap,dtime,freq_cosp, &
8                        ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
[2794]9                        ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
[1925]10                        Nptslmdz,Nlevlmdz,lon,lat, presnivs,overlaplmdz,sunlit, &
[1526]11                        ref_liq,ref_ice,fracTerLic,u_wind,v_wind,phis,phi,ph,p,skt,t, &
[1262]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)
[1925]25! fracTerLic,                           !Fraction terre a convertir en masque
[1262]26! u_wind,v_wind,                        !Vents a 10m ???
27! phi,                                  !Geopotentiel
[1925]28! phis,                                 !Geopotentiel sol
[1262]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
[2428]70
71!!! AI rajouter les nouvelles sorties
[1262]72!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
73
[2428]74!! AI rajouter
[2571]75#include "cosp_defs.h"
[1262]76  USE MOD_COSP_CONSTANTS
77  USE MOD_COSP_TYPES
78  USE MOD_COSP
79  USE mod_phys_lmdz_para
[1327]80  USE mod_grid_phy_lmdz
[1262]81  use ioipsl
82  use iophy
[1925]83  use cosp_output_mod
84  use cosp_output_write_mod
[2428]85!  use MOD_COSP_Modis_Simulator, only : cosp_modis
86
[1262]87  IMPLICIT NONE
88
89  ! Local variables
[1327]90  character(len=64),PARAMETER  :: cosp_input_nl='cosp_input_nl.txt'
91  character(len=64),PARAMETER  :: cosp_output_nl='cosp_output_nl.txt'
[1262]92  integer, save :: isccp_topheight,isccp_topheight_direction,overlap
93  integer,save  :: Ncolumns     ! Number of subcolumns in SCOPS
94  integer, save :: Npoints      ! Number of gridpoints
[1327]95!$OMP THREADPRIVATE(Npoints)
[1262]96  integer, save :: Nlevels      ! Number of levels
97  Integer :: Nptslmdz,Nlevlmdz ! Nb de points issus de physiq.F
98  integer, save :: Nlr          ! Number of levels in statistical outputs
99  integer, save :: Npoints_it   ! Max number of gridpoints to be processed in one iteration
100  integer :: i
101  type(cosp_config),save :: cfg   ! Configuration options
[1327]102!$OMP THREADPRIVATE(cfg)
[1262]103  type(cosp_gridbox) :: gbx ! Gridbox information. Input for COSP
104  type(cosp_subgrid) :: sgx     ! Subgrid outputs
105  type(cosp_sgradar) :: sgradar ! Output from radar simulator
106  type(cosp_sglidar) :: sglidar ! Output from lidar simulator
107  type(cosp_isccp)   :: isccp   ! Output from ISCCP simulator
[2428]108!! AI rajout modis
109  type(cosp_modis)   :: modis   ! Output from MODIS simulator
110!!
[1262]111  type(cosp_misr)    :: misr    ! Output from MISR simulator
[2428]112!! AI rajout rttovs
113!  type(cosp_rttov)   :: rttov   ! Output from RTTOV
114!!
[1262]115  type(cosp_vgrid)   :: vgrid   ! Information on vertical grid of stats
116  type(cosp_radarstats) :: stradar ! Summary statistics from radar simulator
117  type(cosp_lidarstats) :: stlidar ! Summary statistics from lidar simulator
118
119  integer :: t0,t1,count_rate,count_max
[2428]120  integer :: Nlon,Nlat
[1262]121  real,save :: radar_freq,k2,ZenAng,co2,ch4,n2o,co,emsfc_lw
[1327]122!$OMP THREADPRIVATE(emsfc_lw)
[1262]123  integer,dimension(RTTOV_MAX_CHANNELS),save :: Channels
124  real,dimension(RTTOV_MAX_CHANNELS),save :: Surfem
125  integer, save :: surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay
126  integer, save :: Nprmts_max_hydro,Naero,Nprmts_max_aero,lidar_ice_type
127  integer, save :: platform,satellite,Instrument,Nchannels
128  logical, save :: use_vgrid,csat_vgrid,use_precipitation_fluxes,use_reff
129
130! Declaration necessaires pour les sorties IOIPSL
[1925]131  integer :: ii
[2794]132  real    :: ecrit_day,ecrit_hf,ecrit_mth, missing_val
[2137]133  logical :: ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, ok_all_xml
[1925]134
[1262]135  logical, save :: debut_cosp=.true.
[1327]136!$OMP THREADPRIVATE(debut_cosp)
[1526]137
[1262]138!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Input variables from LMDZ-GCM
139  integer                         :: overlaplmdz   !  overlap type: 1=max, 2=rand, 3=max/rand ! cosp input (output lmdz)
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, &
[1526]142                                     zlev,zlev_half,mr_ozone,radliq,radice,dtau_s,dem_s,ref_liq,ref_ice
[1262]143  real,dimension(Nptslmdz,Nlevlmdz) ::  fl_lsrainI,fl_lssnowI,fl_ccrainI,fl_ccsnowI
[1925]144  real,dimension(Nptslmdz)        :: lon,lat,skt,fracTerLic,u_wind,v_wind,phis,sunlit         
[1262]145  real,dimension(Nlevlmdz)        :: presnivs
146  integer                         :: itap,k,ip
147  real                            :: dtime,freq_cosp
[2428]148  real,dimension(2)               :: time_bnds
[1925]149 
[2080]150   namelist/COSP_INPUT/overlap,isccp_topheight,isccp_topheight_direction, &
151              npoints_it,ncolumns,use_vgrid,nlr,csat_vgrid, &
[1262]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
[1327]166  NPoints=Nptslmdz
[1526]167  Nlevels=Nlevlmdz
168 
[1262]169! Lecture du namelist input
[1327]170  CALL read_cosp_input
171
[1262]172! Clefs Outputs
173  call read_cosp_output_nl(cosp_output_nl,cfg)
174
175    if (overlaplmdz.ne.overlap) then
176       print*,'Attention overlaplmdz different de overlap lu dans namelist '
177    endif
178   print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
179
[1925]180  print*,' Cles des differents simulateurs cosp :'
[2428]181  print*,'Lradar_sim,Llidar_sim,Lisccp_sim,Lmisr_sim,Lmodis_sim,Lrttov_sim', &
182          cfg%Lradar_sim,cfg%Llidar_sim,cfg%Lisccp_sim,cfg%Lmisr_sim,cfg%Lmodis_sim,cfg%Lrttov_sim
[1262]183
184  endif ! debut_cosp
185
[2428]186  time_bnds(1) = dtime-dtime/2.
187  time_bnds(2) = dtime+dtime/2.
188
[1925]189!  print*,'Debut phys_cosp itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf ', &
190!          itap,dtime,freq_cosp,ecrit_mth,ecrit_day,ecrit_hf
[1262]191!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
192! Allocate memory for gridbox type
193!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[2835]194! AI mars 2017
195!        print *, 'Allocating memory for gridbox type...'
196
[2428]197!! AI
198!        call construct_cosp_gridbox(dble(itap),radar_freq,surface_radar,use_mie_tables,use_gas_abs,do_ray,melt_lay,k2, &
199!                                    Npoints,Nlevels,Ncolumns,N_HYDRO,Nprmts_max_hydro,Naero,Nprmts_max_aero,Npoints_it, &
200!                                    lidar_ice_type,isccp_topheight,isccp_topheight_direction,overlap,emsfc_lw, &
201!                                    use_precipitation_fluxes,use_reff, &
202!                                    Platform,Satellite,Instrument,Nchannels,ZenAng, &
203!                                    channels(1:Nchannels),surfem(1:Nchannels),co2,ch4,n2o,co,gbx)
204! Surafce emissivity
205        emsfc_lw = 1.
[1262]206
[2428]207        call construct_cosp_gridbox(dtime,time_bnds,radar_freq,surface_radar,use_mie_tables,use_gas_abs, &
208                                    do_ray,melt_lay,k2, &
[1262]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
[2835]219!        print *, 'Populating input structure...'
[1262]220        gbx%longitude = lon
221        gbx%latitude = lat
222
223        gbx%p = p !
224        gbx%ph = ph
[1526]225        gbx%zlev = phi/9.81
[1262]226
[1526]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
[1262]235
236        gbx%T = T
237        gbx%q = rh*100.
238        gbx%sh = sh
[1925]239! On ne veut pas que cosp distingue les nuages stratiformes et convectifs
240! on passe les contenus totaux (conv+strat)
241        gbx%cca = 0. !convective_cloud_amount (1)
[1262]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
[1925]258! sunlit calcule a partir de la fraction d ensoleillement par jour
[2428]259!      do ip = 1, Npoints
260!        if (sunlit(ip).le.0.) then
261!           gbx%sunlit(ip)=0.
262!        else
263!           gbx%sunlit(ip)=1.
264!        endif
265!      enddo
266       gbx%sunlit=sunlit
[1925]267
[1262]268! A voir l equivalent LMDZ
269  mr_ccliq = 0.0
270  mr_ccice = 0.0
271        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
272        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
273        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
274        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
275! A revoir
276        fl_lsrain = fl_lsrainI + fl_ccrainI
277        fl_lssnow = fl_lssnowI + fl_ccsnowI
278        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
279        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
280!  A voir l equivalent LMDZ
281        fl_lsgrpl=0.
282        fl_ccsnow = 0.
283        fl_ccrain = 0.
284        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
285        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
286        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
287
288     gbx%Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
289     gbx%Reff(:,:,I_LSCICE) = ref_ice*1e-6
[2428]290!! AI A revoir
[1262]291     gbx%Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
292     gbx%Reff(:,:,I_CVCICE) = ref_ice*1e-6
293
294        ! ISCCP simulator
295        gbx%dtau_s   = dtau_s
296        gbx%dtau_c   = 0.
297        gbx%dem_s    = dem_s
298        gbx%dem_c    = 0.
299
300!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
301        ! Define new vertical grid
302!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[2835]303!        print *, 'Defining new vertical grid...'
[1262]304        call construct_cosp_vgrid(gbx,Nlr,use_vgrid,csat_vgrid,vgrid)
305
306!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
307       ! Allocate memory for other types
308!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[2835]309!        print *, 'Allocating memory for other types...'
[1262]310        call construct_cosp_subgrid(Npoints, Ncolumns, Nlevels, sgx)
311        call construct_cosp_sgradar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,sgradar)
312        call construct_cosp_radarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,stradar)
313        call construct_cosp_sglidar(cfg,Npoints,Ncolumns,Nlevels,N_HYDRO,PARASOL_NREFL,sglidar)
314        call construct_cosp_lidarstats(cfg,Npoints,Ncolumns,vgrid%Nlvgrid,N_HYDRO,PARASOL_NREFL,stlidar)
315        call construct_cosp_isccp(cfg,Npoints,Ncolumns,Nlevels,isccp)
[2428]316!! AI rajout
317        call construct_cosp_modis(cfg,Npoints,modis)
318!!
[1262]319        call construct_cosp_misr(cfg,Npoints,misr)
[2428]320!        call construct_cosp_rttov(cfg,Npoints,Nchannels,rttov)
[1925]321
322!+++++++++++++ Open output files and define output files axis !+++++++++++++
323     if (debut_cosp) then
324
[2080]325      !$OMP MASTER
[1925]326        print *, ' Open outpts files and define axis'
327        call cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp, &
[2137]328                              ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml, &
[2822]329                              ecrit_mth, ecrit_day, ecrit_hf, use_vgrid, vgrid, stlidar)
[2080]330      !$OMP END MASTER
331      !$OMP BARRIER
[1925]332        debut_cosp=.false.
333      endif ! debut_cosp
[1262]334!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
335        ! Call simulator
336!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[2835]337!        print *, 'Calling simulator...'
[2428]338!! AI
339!        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,stradar,stlidar)
340!#ifdef RTTOV
341!        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,modis,rttov,stradar,stlidar)
342!#else
343        call cosp(overlap,Ncolumns,cfg,vgrid,gbx,sgx,sgradar,sglidar,isccp,misr,modis,stradar,stlidar)
344!#endif
345!!
346
[1262]347!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
348
[1925]349!!!!!!!!!!!!!!!!!! Ecreture des sorties Cosp !!!!!!!!!!!!!!r!!!!!!:!!!!!
[2428]350
[2835]351!       print *, 'Calling write output'
[2794]352        call cosp_output_write(Nlevlmdz, Npoints, Ncolumns, itap, dtime, freq_COSP, missing_val, &
[2428]353                               cfg, gbx, vgrid, sglidar, sgradar, stlidar, stradar, &
354                               isccp, misr, modis)
[1262]355
356!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
357        ! Deallocate memory in derived types
358!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[2835]359!        print *, 'Deallocating memory...'
[1262]360        call free_cosp_gridbox(gbx)
361        call free_cosp_subgrid(sgx)
362        call free_cosp_sgradar(sgradar)
363        call free_cosp_radarstats(stradar)
364        call free_cosp_sglidar(sglidar)
365        call free_cosp_lidarstats(stlidar)
366        call free_cosp_isccp(isccp)
367        call free_cosp_misr(misr)
[2428]368!! AI
369        call free_cosp_modis(modis)
370!        call free_cosp_rttov(rttov)
371!!
[1262]372        call free_cosp_vgrid(vgrid) 
373 
374!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
375  ! Time in s. Only for testing purposes
376!  call system_clock(t1,count_rate,count_max)
377!  print *,(t1-t0)*1.0/count_rate
[1327]378 
379  CONTAINS
380 
381  SUBROUTINE read_cosp_input
[1262]382   
[1327]383    IF (is_master) THEN
384      OPEN(10,file=cosp_input_nl,status='old')
385      READ(10,nml=cosp_input)
386      CLOSE(10)
387    ENDIF
388    CALL bcast(overlap)
389    CALL bcast(isccp_topheight)
390    CALL bcast(isccp_topheight_direction)
391    CALL bcast(npoints_it)
392    CALL bcast(ncolumns)
393    CALL bcast(use_vgrid)
394    CALL bcast(nlr)
395    CALL bcast(csat_vgrid)
396    CALL bcast(radar_freq)
397    CALL bcast(surface_radar)
398    CALL bcast(use_mie_tables)
399    CALL bcast(use_gas_abs)
400    CALL bcast(do_ray)
401    CALL bcast(melt_lay)
402    CALL bcast(k2)
403    CALL bcast(Nprmts_max_hydro)
404    CALL bcast(Naero)
405    CALL bcast(Nprmts_max_aero)
406    CALL bcast(lidar_ice_type)
407    CALL bcast(use_precipitation_fluxes)
408    CALL bcast(use_reff)
409    CALL bcast(platform)
410    CALL bcast(satellite)
411    CALL bcast(Instrument)
412    CALL bcast(Nchannels)
413    CALL bcast(Channels)
414    CALL bcast(Surfem)
415    CALL bcast(ZenAng)
416    CALL bcast(co2)
417    CALL bcast(ch4)
418    CALL bcast(n2o)
419    CALL bcast(co)
420!$OMP BARRIER 
421  END SUBROUTINE read_cosp_input
422
[1262]423end subroutine phys_cosp
Note: See TracBrowser for help on using the repository browser.