source: LMDZ6/trunk/libf/phylmd/cosp/phys_cosp.F90 @ 4640

Last change on this file since 4640 was 4619, checked in by yann meurdesoif, 13 months ago

Suppress usage of preprocessing key CPP_XIOS.
Wrapper file is used to suppress XIOS symbol when xios is not linked and not used (-io ioipsl)
The CPP_XIOS key is replaced in model by "using_xios" boolean variable to switch between IOIPSL or XIOS output.

YM

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