source: LMDZ6/branches/LMDZ_ECRad/libf/phylmd/ecrad/driver/ecrad_ifs_driver.F90 @ 4728

Last change on this file since 4728 was 4728, checked in by idelkadi, 10 months ago

Update of ecrad in the LMDZ_ECRad branch of LMDZ:

  • version 1.6.1 of ecrad
  • files are no longer grouped in the same ecrad directory.
  • the structure of ecrad offline is preserved to facilitate updating in LMDZ
  • cfg.bld modified to take into account the new added subdirectories.
  • the interface routines and those added in ecrad are moved to the phylmd directory
File size: 19.4 KB
Line 
1! ecrad_ifs_driver.F90 - Driver for offline ECRAD radiation scheme
2!
3! (C) Copyright 2014- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15! ECRAD is the radiation scheme used in the ECMWF Integrated
16! Forecasting System in cycle 43R3 and later. Several solvers are
17! available, including McICA, Tripleclouds and SPARTACUS (the Speedy
18! Algorithm for Radiative Transfer through Cloud Sides, a modification
19! of the two-stream formulation of shortwave and longwave radiative
20! transfer to account for 3D radiative effects). Gas optical
21! properties are provided by the RRTM-G gas optics scheme.
22
23! This program takes three arguments:
24! 1) Namelist file to configure the radiation calculation, but note
25!    that only the radiation_config group is read
26! 2) Name of a NetCDF file containing one or more atmospheric profiles
27! 3) Name of output NetCDF file
28!
29! This version uses the infrastructure of the IFS, such as computing
30! effective radius and cloud overlap from latitude and other
31! variables. To configure ecRad in this version you need to edit
32! ifs/yoerad.F90 in the ecRad package, but these options can be
33! overridden with the "radiation" namelist. This file requires the
34! input data to have compatible settings, e.g. the right number of
35! aerosol variables, and surface albedo/emissivity bands; a test file
36! satisfying this requirement is test/ifs/ecrad_meridian.nc in the
37! ecRad package.
38!
39! Note that the purpose of this file is simply to demonstrate the use
40! of the setup_radiation_scheme and radiation_scheme routines; all the
41! rest is using the offline ecRad driver containers to read a NetCDF
42! file to memory and pass it into these routines.
43
44program ecrad_ifs_driver
45
46  ! --------------------------------------------------------
47  ! Section 1: Declarations
48  ! --------------------------------------------------------
49  use parkind1,                 only : jprb, jprd ! Working/double precision
50
51  use radiation_io,             only : nulout
52  use radiation_single_level,   only : single_level_type
53  use radiation_thermodynamics, only : thermodynamics_type
54  use radiation_gas,            only : gas_type, IMassMixingRatio, &
55       &   IH2O, ICO2, IO3, IN2O, INO2, ICO, ICH4, IO2, ICFC11, ICFC12, &
56       &   IHCFC22, ICCl4
57  use radiation_cloud,          only : cloud_type
58  use radiation_aerosol,        only : aerosol_type
59  use radiation_flux,           only : flux_type
60  use radiation_save,           only : save_net_fluxes
61  use radiation_setup,          only : tradiation, setup_radiation_scheme
62  use radiation_constants,      only : Pi
63  use ecrad_driver_config,      only : driver_config_type
64  use ecrad_driver_read_input,  only : read_input
65  use easy_netcdf
66
67  implicit none
68
69#include "radiation_scheme.intfb.h"
70
71  ! The NetCDF file containing the input profiles
72  type(netcdf_file)         :: file
73
74  ! Configuration for the radiation scheme, IFS style
75  type(tradiation)          :: yradiation
76
77  ! Derived types for the inputs to the radiation scheme
78  type(single_level_type)   :: single_level
79  type(thermodynamics_type) :: thermodynamics
80  type(gas_type)            :: gas
81  type(cloud_type)          :: cloud
82  type(aerosol_type)        :: aerosol
83
84  ! Configuration specific to this driver
85  type(driver_config_type)  :: driver_config
86
87  ! Derived type containing outputs from the radiation scheme
88  type(flux_type)           :: flux
89
90  ! Additional arrays passed to radiation_scheme
91  real(jprb), allocatable, dimension(:) :: ccn_land, ccn_sea, sin_latitude, longitude_rad, land_frac
92  real(jprb), allocatable, dimension(:,:) :: pressure_fl, temperature_fl, zeros
93  real(jprb), allocatable, dimension(:,:,:) :: tegen_aerosol
94  real(jprb), allocatable, dimension(:) :: flux_sw_direct_normal, flux_uv, flux_par, flux_par_clear, &
95       &  flux_incoming, emissivity_out
96  real(jprb), allocatable, dimension(:,:) :: flux_diffuse_band, flux_direct_band
97  real(jprb), allocatable, dimension(:,:) :: cloud_fraction, cloud_q_liq, cloud_q_ice
98
99  integer :: ncol, nlev         ! Number of columns and levels
100  integer :: istartcol, iendcol ! Range of columns to process
101
102  ! Name of file names specified on command line
103  character(len=512) :: file_name
104  integer            :: istatus ! Result of command_argument_count
105
106  ! For parallel processing of multiple blocks
107  integer :: jblock, nblock ! Block loop index and number
108
109#ifndef NO_OPENMP
110  ! OpenMP functions
111  integer, external :: omp_get_thread_num
112  real(kind=jprd), external :: omp_get_wtime
113  ! Start/stop time in seconds
114  real(kind=jprd) :: tstart, tstop
115#endif
116
117  ! For demonstration of get_sw_weights later on
118  ! Ultraviolet weightings
119  !integer    :: nweight_uv
120  !integer    :: iband_uv(100)
121  !real(jprb) :: weight_uv(100)
122  ! Photosynthetically active radiation weightings
123  !integer    :: nweight_par
124  !integer    :: iband_par(100)
125  !real(jprb) :: weight_par(100)
126
127  ! Loop index for repeats (for benchmarking)
128  integer :: jrepeat
129
130  ! Are any variables out of bounds?
131  logical :: is_out_of_bounds
132
133!  integer    :: iband(20), nweights
134!  real(jprb) :: weight(20)
135
136
137  ! --------------------------------------------------------
138  ! Section 2: Configure
139  ! --------------------------------------------------------
140
141  ! Check program called with correct number of arguments
142  if (command_argument_count() < 3) then
143    stop 'Usage: ecrad config.nam input_file.nc output_file.nc'
144  end if
145
146  ! Use namelist to configure the radiation calculation
147  call get_command_argument(1, file_name, status=istatus)
148  if (istatus /= 0) then
149    stop 'Failed to read name of namelist file as string of length < 512'
150  end if
151
152  ! Read "radiation_driver" namelist into radiation driver config type
153  call driver_config%read(file_name)
154
155  if (driver_config%iverbose >= 2) then
156    write(nulout,'(a)') '-------------------------- OFFLINE ECRAD RADIATION SCHEME --------------------------'
157    write(nulout,'(a)') 'Copyright (C) 2014- ECMWF'
158    write(nulout,'(a)') 'Contact: Robin Hogan (r.j.hogan@ecmwf.int)'
159#ifdef PARKIND1_SINGLE
160    write(nulout,'(a)') 'Floating-point precision: single'
161#else
162    write(nulout,'(a)') 'Floating-point precision: double'
163#endif
164  end if
165
166  ! Albedo/emissivity intervals may be specified like this
167  !call config%define_sw_albedo_intervals(6, &
168  !     &  [0.25e-6_jprb, 0.44e-6_jprb, 0.69e-6_jprb, &
169  !     &     1.19_jprb, 2.38e-6_jprb], [1,2,3,4,5,6], &
170  !     &   do_nearest=.false.)
171  !call config%define_lw_emiss_intervals(3, &
172  !     &  [8.0e-6_jprb, 13.0e-6_jprb], [1,2,1], &
173  !     &   do_nearest=.false.)
174
175  ! If monochromatic aerosol properties are required, then the
176  ! wavelengths can be specified (in metres) as follows - these can be
177  ! whatever you like for the general aerosol optics, but must match
178  ! the monochromatic values in the aerosol input file for the older
179  ! aerosol optics
180  !call config%set_aerosol_wavelength_mono( &
181  !     &  [3.4e-07_jprb, 3.55e-07_jprb, 3.8e-07_jprb, 4.0e-07_jprb, 4.4e-07_jprb, &
182  !     &   4.69e-07_jprb, 5.0e-07_jprb, 5.32e-07_jprb, 5.5e-07_jprb, 6.45e-07_jprb, &
183  !     &   6.7e-07_jprb, 8.0e-07_jprb, 8.58e-07_jprb, 8.65e-07_jprb, 1.02e-06_jprb, &
184  !     &   1.064e-06_jprb, 1.24e-06_jprb, 1.64e-06_jprb, 2.13e-06_jprb, 1.0e-05_jprb])
185
186  call yradiation%rad_config%read(file_name=file_name)
187
188  ! Setup aerosols
189  if (yradiation%rad_config%use_aerosols) then
190    yradiation%yrerad%naermacc = 1 ! MACC-derived aerosol climatology on a NMCLAT x NMCLON grid
191  else
192    yradiation%yrerad%naermacc = 0
193  endif
194
195  ! Setup the radiation scheme: load the coefficients for gas and
196  ! cloud optics, currently from RRTMG
197  call setup_radiation_scheme(yradiation, .true., file_name=file_name)
198  ! Or call without specifying the namelist filename, in which case
199  ! the default settings are from yoerad.F90
200  !call setup_radiation_scheme(yradiation, .true.)
201
202  ! Demonstration of how to get weights for UV and PAR fluxes
203  !if (config%do_sw) then
204  !  call config%get_sw_weights(0.2e-6_jprb, 0.4415e-6_jprb,&
205  !       &  nweight_uv, iband_uv, weight_uv,&
206  !       &  'ultraviolet')
207  !  call config%get_sw_weights(0.4e-6_jprb, 0.7e-6_jprb,&
208  !       &  nweight_par, iband_par, weight_par,&
209  !       &  'photosynthetically active radiation, PAR')
210  !end if
211
212  ! --------------------------------------------------------
213  ! Section 3: Read input data file
214  ! --------------------------------------------------------
215
216  ! Get NetCDF input file name
217  call get_command_argument(2, file_name, status=istatus)
218  if (istatus /= 0) then
219    stop 'Failed to read name of input NetCDF file as string of length < 512'
220  end if
221
222  ! Open the file and configure the way it is read
223  call file%open(trim(file_name), iverbose=driver_config%iverbose)
224
225  ! Get NetCDF output file name
226  call get_command_argument(3, file_name, status=istatus)
227  if (istatus /= 0) then
228    stop 'Failed to read name of output NetCDF file as string of length < 512'
229  end if
230
231  ! 2D arrays are assumed to be stored in the file with height varying
232  ! more rapidly than column index. Specifying "true" here transposes
233  ! all 2D arrays so that the column index varies fastest within the
234  ! program.
235  call file%transpose_matrices(.true.)
236
237  ! Read input variables from NetCDF file, noting that cloud overlap
238  ! and effective radius are ignored
239  call read_input(file, yradiation%rad_config, driver_config, ncol, nlev, &
240       &          single_level, thermodynamics, &
241       &          gas, cloud, aerosol)
242
243  ! Latitude is used for cloud overlap and ice effective radius
244  if (file%exists('lat')) then
245    call file%get('lat', sin_latitude)
246    sin_latitude = sin(sin_latitude * Pi/180.0_jprb)
247  else
248    allocate(sin_latitude(ncol))
249    sin_latitude = 0.0_jprb
250  end if
251
252  if (file%exists('lon')) then
253    call file%get('lon', longitude_rad)
254    longitude_rad = longitude_rad * Pi/180.0_jprb
255  else
256    allocate(longitude_rad(ncol))
257    longitude_rad = 0.0_jprb
258  end if
259
260  ! Close input file
261  call file%close()
262
263  ! Convert gas units to mass-mixing ratio
264  call gas%set_units(IMassMixingRatio)
265
266  ! Compute seed from skin temperature residual
267  !  single_level%iseed = int(1.0e9*(single_level%skin_temperature &
268  !       &                            -int(single_level%skin_temperature)))
269
270  ! Set first and last columns to process
271  if (driver_config%iendcol < 1 .or. driver_config%iendcol > ncol) then
272    driver_config%iendcol = ncol
273  end if
274
275  if (driver_config%istartcol > driver_config%iendcol) then
276    write(nulout,'(a,i0,a,i0,a,i0,a)') '*** Error: requested column range (', &
277         &  driver_config%istartcol, &
278         &  ' to ', driver_config%iendcol, ') is out of the range in the data (1 to ', &
279         &  ncol, ')'
280    stop 1
281  end if
282
283  ! --------------------------------------------------------
284  ! Section 4: Call radiation scheme
285  ! --------------------------------------------------------
286
287  ! Compute saturation with respect to liquid (needed for aerosol
288  ! hydration) call
289  !  call thermodynamics%calc_saturation_wrt_liquid(driver_config%istartcol,driver_config%iendcol)
290
291  ! Check inputs are within physical bounds, printing message if not
292  is_out_of_bounds =     gas%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol, &
293       &                                            driver_config%do_correct_unphysical_inputs) &
294       & .or.   single_level%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol, &
295       &                                            driver_config%do_correct_unphysical_inputs) &
296       & .or. thermodynamics%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol, &
297       &                                            driver_config%do_correct_unphysical_inputs) &
298       & .or.          cloud%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol, &
299       &                                            driver_config%do_correct_unphysical_inputs) &
300       & .or.        aerosol%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol, &
301       &                                            driver_config%do_correct_unphysical_inputs)
302
303  ! Allocate memory for the flux profiles, which may include arrays
304  ! of dimension n_bands_sw/n_bands_lw, so must be called after
305  ! setup_radiation
306  call flux%allocate(yradiation%rad_config, 1, ncol, nlev)
307
308  ! set relevant fluxes to zero
309  flux%lw_up(:,:) = 0._jprb
310  flux%lw_dn(:,:) = 0._jprb
311  flux%sw_up(:,:) = 0._jprb
312  flux%sw_dn(:,:) = 0._jprb
313  flux%sw_dn_direct(:,:) = 0._jprb
314  flux%lw_up_clear(:,:) = 0._jprb
315  flux%lw_dn_clear(:,:) = 0._jprb
316  flux%sw_up_clear(:,:) = 0._jprb
317  flux%sw_dn_clear(:,:) = 0._jprb
318  flux%sw_dn_direct_clear(:,:) = 0._jprb
319
320  flux%lw_dn_surf_canopy(:,:) = 0._jprb
321  flux%sw_dn_diffuse_surf_canopy(:,:) = 0._jprb
322  flux%sw_dn_direct_surf_canopy(:,:) = 0._jprb
323  flux%lw_derivatives(:,:) = 0._jprb
324
325  ! Allocate memory for additional arrays
326  allocate(ccn_land(ncol))
327  allocate(ccn_sea(ncol))
328  allocate(land_frac(ncol))
329  allocate(pressure_fl(ncol,nlev))
330  allocate(temperature_fl(ncol,nlev))
331  allocate(zeros(ncol,nlev))
332  allocate(tegen_aerosol(ncol,6,nlev))
333  allocate(flux_sw_direct_normal(ncol))
334  allocate(flux_uv(ncol))
335  allocate(flux_par(ncol))
336  allocate(flux_par_clear(ncol))
337  allocate(flux_incoming(ncol))
338  allocate(emissivity_out(ncol))
339  allocate(flux_diffuse_band(ncol,yradiation%yrerad%nsw))
340  allocate(flux_direct_band(ncol,yradiation%yrerad%nsw))
341  allocate(cloud_fraction(ncol,nlev))
342  allocate(cloud_q_liq(ncol,nlev))
343  allocate(cloud_q_ice(ncol,nlev))
344
345  ccn_land = yradiation%yrerad%rccnlnd
346  ccn_sea = yradiation%yrerad%rccnsea
347  tegen_aerosol = 0.0_jprb
348  pressure_fl = 0.5_jprb * (thermodynamics%pressure_hl(:,1:nlev)+thermodynamics%pressure_hl(:,2:nlev+1))
349  temperature_fl = 0.5_jprb * (thermodynamics%temperature_hl(:,1:nlev)+thermodynamics%temperature_hl(:,2:nlev+1))
350  zeros = 0.0_jprb ! Dummy snow/rain water mixing ratios
351
352  if (yradiation%rad_config%do_clouds) then
353    cloud_fraction = cloud%fraction
354    cloud_q_liq = cloud%q_liq
355    cloud_q_ice = cloud%q_ice
356  else
357    cloud_fraction = 0.0_jprb
358    cloud_q_liq = 0.0_jprb
359    cloud_q_ice = 0.0_jprb
360  endif
361
362  if (driver_config%iverbose >= 2) then
363    write(nulout,'(a)')  'Performing radiative transfer calculations'
364  end if
365
366  ! Option of repeating calculation multiple time for more accurate
367  ! profiling
368#ifndef NO_OPENMP
369  tstart = omp_get_wtime()
370#endif
371  do jrepeat = 1,driver_config%nrepeat
372
373!    if (driver_config%do_parallel) then
374      ! Run radiation scheme over blocks of columns in parallel
375
376      ! Compute number of blocks to process
377      nblock = (driver_config%iendcol - driver_config%istartcol &
378           &  + driver_config%nblocksize) / driver_config%nblocksize
379
380      !$OMP PARALLEL DO PRIVATE(istartcol, iendcol) SCHEDULE(RUNTIME)
381      do jblock = 1, nblock
382        ! Specify the range of columns to process.
383        istartcol = (jblock-1) * driver_config%nblocksize &
384             &    + driver_config%istartcol
385        iendcol = min(istartcol + driver_config%nblocksize - 1, &
386             &        driver_config%iendcol)
387
388        if (driver_config%iverbose >= 3) then
389#ifndef NO_OPENMP
390          write(nulout,'(a,i0,a,i0,a,i0)')  'Thread ', omp_get_thread_num(), &
391               &  ' processing columns ', istartcol, '-', iendcol
392#else
393          write(nulout,'(a,i0,a,i0)')  'Processing columns ', istartcol, '-', iendcol
394#endif
395        end if
396
397        ! Call the ECRAD radiation scheme; note that we are simply
398        ! passing arrays in rather than ecRad structures, which are
399        ! used here just for convenience
400        call radiation_scheme(yradiation, istartcol, iendcol, ncol, nlev, size(aerosol%mixing_ratio,3), &
401             &  single_level%solar_irradiance, single_level%cos_sza, single_level%skin_temperature, &
402             &  single_level%sw_albedo, single_level%sw_albedo_direct, single_level%lw_emissivity, &
403             &  ccn_land, ccn_sea, longitude_rad, sin_latitude, land_frac, pressure_fl, temperature_fl, &
404             &  thermodynamics%pressure_hl, thermodynamics%temperature_hl, &
405             &  gas%mixing_ratio(:,:,IH2O), gas%mixing_ratio(:,:,ICO2), &
406             &  gas%mixing_ratio(:,:,ICH4), gas%mixing_ratio(:,:,IN2O), gas%mixing_ratio(:,:,INO2), &
407             &  gas%mixing_ratio(:,:,ICFC11), gas%mixing_ratio(:,:,ICFC12), gas%mixing_ratio(:,:,IHCFC22), &
408             &  gas%mixing_ratio(:,:,ICCl4), gas%mixing_ratio(:,:,IO3), cloud_fraction, cloud_q_liq, &
409             &  cloud_q_ice, zeros, zeros, tegen_aerosol, aerosol%mixing_ratio, flux%sw_up, flux%lw_up, &
410             &  flux%sw_up_clear, flux%lw_up_clear, flux%sw_dn(:,nlev+1), flux%lw_dn(:,nlev+1), &
411             &  flux%sw_dn_clear(:,nlev+1), flux%lw_dn_clear(:,nlev+1), &
412             &  flux%sw_dn_direct(:,nlev+1), flux%sw_dn_direct_clear(:,nlev+1), flux_sw_direct_normal, &
413             &  flux_uv, flux_par, &
414             &  flux_par_clear, flux%sw_dn(:,1), emissivity_out, flux%lw_derivatives, flux_diffuse_band, &
415             &  flux_direct_band)
416      end do
417      !$OMP END PARALLEL DO
418
419!    else
420      ! Run radiation scheme serially
421!      if (driver_config%iverbose >= 3) then
422!        write(nulout,'(a,i0,a)')  'Processing ', ncol, ' columns'
423!      end if
424
425      ! Call the ECRAD radiation scheme
426!      call radiation_scheme(ncol, nlev, driver_config%istartcol, driver_config%iendcol, &
427!           &  config, single_level, thermodynamics, gas, cloud, aerosol, flux)
428
429!    end if
430
431  end do
432
433  ! "up" fluxes are actually net fluxes at this point - we modify the
434  ! upwelling flux so that net=dn-up, while the TOA and surface
435  ! downwelling fluxes are correct.
436  flux%sw_up = -flux%sw_up
437  flux%sw_up(:,1) = flux%sw_up(:,1)+flux%sw_dn(:,1)
438  flux%sw_up(:,nlev+1) = flux%sw_up(:,nlev+1)+flux%sw_dn(:,nlev+1)
439
440  flux%lw_up = -flux%lw_up
441  flux%lw_up(:,1) = flux%lw_up(:,1)+flux%lw_dn(:,1)
442  flux%lw_up(:,nlev+1) = flux%lw_up(:,nlev+1)+flux%lw_dn(:,nlev+1)
443
444  flux%sw_up_clear = -flux%sw_up_clear
445  flux%sw_up_clear(:,1) = flux%sw_up_clear(:,1)+flux%sw_dn_clear(:,1)
446  flux%sw_up_clear(:,nlev+1) = flux%sw_up_clear(:,nlev+1)+flux%sw_dn_clear(:,nlev+1)
447
448  flux%lw_up_clear = -flux%lw_up_clear
449  flux%lw_up_clear(:,1) = flux%lw_up_clear(:,1)+flux%lw_dn_clear(:,1)
450  flux%lw_up_clear(:,nlev+1) = flux%lw_up_clear(:,nlev+1)+flux%lw_dn_clear(:,nlev+1)
451
452#ifndef NO_OPENMP
453  tstop = omp_get_wtime()
454  write(nulout, '(a,g12.5,a)') 'Time elapsed in radiative transfer: ', tstop-tstart, ' seconds'
455#endif
456
457  ! --------------------------------------------------------
458  ! Section 5: Check and save output
459  ! --------------------------------------------------------
460
461  ! This is unreliable because only the net fluxes are valid:
462  !is_out_of_bounds = flux%out_of_physical_bounds(driver_config%istartcol, driver_config%iendcol)
463
464  ! Store the fluxes in the output file
465  yradiation%rad_config%do_surface_sw_spectral_flux = .false.
466  yradiation%rad_config%do_canopy_fluxes_sw = .false.
467  yradiation%rad_config%do_canopy_fluxes_lw = .false.
468
469  call save_net_fluxes(file_name, yradiation%rad_config, thermodynamics, flux, &
470       &   iverbose=driver_config%iverbose, is_hdf5_file=driver_config%do_write_hdf5, &
471       &   experiment_name=driver_config%experiment_name, &
472       &   is_double_precision=driver_config%do_write_double_precision)
473
474  if (driver_config%iverbose >= 2) then
475    write(nulout,'(a)') '------------------------------------------------------------------------------------'
476  end if
477
478end program ecrad_ifs_driver
Note: See TracBrowser for help on using the repository browser.