source: palm/trunk/SOURCE/radiation_model_mod.f90 @ 1826

Last change on this file since 1826 was 1826, checked in by maronga, 8 years ago

further modularization of radiation model and plant canopy model

  • Property svn:keywords set to Id
File size: 87.4 KB
RevLine 
[1826]1!> @file radiation_model_mod.f90
[1496]2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
[1818]16! Copyright 1997-2016 Leibniz Universitaet Hannover
[1496]17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
[1826]21! Further modularization.
[1789]22!
[1496]23! Former revisions:
24! -----------------
25! $Id: radiation_model_mod.f90 1826 2016-04-07 12:01:39Z maronga $
26!
[1789]27! 1788 2016-03-10 11:01:04Z maronga
28! Added new albedo class for pavements / roads.
29!
[1784]30! 1783 2016-03-06 18:36:17Z raasch
31! palm-netcdf-module removed in order to avoid a circular module dependency,
32! netcdf-variables moved to netcdf-module, new routine netcdf_handle_error_rad
33! added
34!
[1758]35! 1757 2016-02-22 15:49:32Z maronga
36! Added parameter unscheduled_radiation_calls. Bugfix: interpolation of sounding
37! profiles for pressure and temperature above the LES domain.
38!
[1710]39! 1709 2015-11-04 14:47:01Z maronga
40! Bugfix: set initial value for rrtm_lwuflx_dt to zero, small formatting
41! corrections
42!
[1702]43! 1701 2015-11-02 07:43:04Z maronga
44! Bugfixes: wrong index for output of timeseries, setting of nz_snd_end
45!
[1692]46! 1691 2015-10-26 16:17:44Z maronga
47! Added option for spin-up runs without radiation (skip_time_do_radiation). Bugfix
48! in calculation of pressure profiles. Bugfix in calculation of trace gas profiles.
49! Added output of radiative heating rates.
50!
[1683]51! 1682 2015-10-07 23:56:08Z knoop
52! Code annotations made doxygen readable
53!
[1607]54! 1606 2015-06-29 10:43:37Z maronga
55! Added preprocessor directive __netcdf to allow for compiling without netCDF.
56! Note, however, that RRTMG cannot be used without netCDF.
57!
[1591]58! 1590 2015-05-08 13:56:27Z maronga
59! Bugfix: definition of character strings requires same length for all elements
60!
[1588]61! 1587 2015-05-04 14:19:01Z maronga
62! Added albedo class for snow
63!
[1586]64! 1585 2015-04-30 07:05:52Z maronga
65! Added support for RRTMG
66!
[1572]67! 1571 2015-03-12 16:12:49Z maronga
68! Added missing KIND attribute. Removed upper-case variable names
69!
[1552]70! 1551 2015-03-03 14:18:16Z maronga
71! Added support for data output. Various variables have been renamed. Added
72! interface for different radiation schemes (currently: clear-sky, constant, and
73! RRTM (not yet implemented).
74!
[1497]75! 1496 2014-12-02 17:25:50Z maronga
76! Initial revision
77!
[1496]78!
79! Description:
80! ------------
[1682]81!> Radiation models and interfaces
[1826]82!> @todo move variable definitions used in radiation_init only to the subroutine
[1682]83!>       as they are no longer required after initialization.
84!> @todo Output of full column vertical profiles used in RRTMG
85!> @todo Output of other rrtm arrays (such as volume mixing ratios)
86!> @todo Adapt for use with topography
87!>
88!> @note Many variables have a leading dummy dimension (0:0) in order to
89!>       match the assume-size shape expected by the RRTMG model.
[1496]90!------------------------------------------------------------------------------!
[1682]91 MODULE radiation_model_mod
92 
[1496]93    USE arrays_3d,                                                             &
[1691]94        ONLY:  dzw, hyp, pt, q, ql, zw
[1496]95
[1585]96    USE cloud_parameters,                                                      &
[1691]97        ONLY:  cp, l_d_cp, nc_const, rho_l, sigma_gc 
[1585]98
99    USE constants,                                                             &
100        ONLY:  pi
101
[1496]102    USE control_parameters,                                                    &
[1585]103        ONLY:  cloud_droplets, cloud_physics, g, initializing_actions,         &
[1691]104               large_scale_forcing, lsf_surf, phi, pt_surface, rho_surface,    &
[1585]105               surface_pressure, time_since_reference_point
[1496]106
107    USE indices,                                                               &
[1585]108        ONLY:  nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb_s_inner, nzb, nzt
[1496]109
110    USE kinds
111
[1606]112#if defined ( __netcdf )
[1783]113    USE NETCDF
[1606]114#endif
[1585]115
116#if defined ( __rrtmg )
[1788]117
[1585]118    USE parrrsw,                                                               &
119        ONLY:  naerec, nbndsw
[1551]120
[1585]121    USE parrrtm,                                                               &
122        ONLY:  nbndlw
123
124    USE rrtmg_lw_init,                                                         &
125        ONLY:  rrtmg_lw_ini
126
127    USE rrtmg_sw_init,                                                         &
128        ONLY:  rrtmg_sw_ini
129
130    USE rrtmg_lw_rad,                                                          &
131        ONLY:  rrtmg_lw
132
133    USE rrtmg_sw_rad,                                                          &
134        ONLY:  rrtmg_sw
135#endif
136
137
138
[1496]139    IMPLICIT NONE
140
[1585]141    CHARACTER(10) :: radiation_scheme = 'clear-sky' ! 'constant', 'clear-sky', or 'rrtmg'
[1551]142
[1585]143!
144!-- Predefined Land surface classes (albedo_type) after Briegleb (1992)
[1788]145    CHARACTER(37), DIMENSION(0:17), PARAMETER :: albedo_type_name = (/      &
[1590]146                                   'user defined                         ', & !  0
147                                   'ocean                                ', & !  1
148                                   'mixed farming, tall grassland        ', & !  2
149                                   'tall/medium grassland                ', & !  3
150                                   'evergreen shrubland                  ', & !  4
151                                   'short grassland/meadow/shrubland     ', & !  5
152                                   'evergreen needleleaf forest          ', & !  6
153                                   'mixed deciduous evergreen forest     ', & !  7
154                                   'deciduous forest                     ', & !  8
155                                   'tropical evergreen broadleaved forest', & !  9
156                                   'medium/tall grassland/woodland       ', & ! 10
157                                   'desert, sandy                        ', & ! 11
158                                   'desert, rocky                        ', & ! 12
159                                   'tundra                               ', & ! 13
160                                   'land ice                             ', & ! 14
161                                   'sea ice                              ', & ! 15
[1788]162                                   'snow                                 ', & ! 16
163                                   'pavement/roads                       '  & ! 17
[1585]164                                                         /)
[1496]165
[1682]166    INTEGER(iwp) :: albedo_type  = 5,    & !< Albedo surface type (default: short grassland)
167                    day,                 & !< current day of the year
168                    day_init     = 172,  & !< day of the year at model start (21/06)
169                    dots_rad     = 0       !< starting index for timeseries output
[1496]170
[1757]171    LOGICAL ::  unscheduled_radiation_calls = .TRUE., & !< flag parameter indicating whether additional calls of the radiation code are allowed
172                constant_albedo = .FALSE.,            & !< flag parameter indicating whether the albedo may change depending on zenith
173                force_radiation_call = .FALSE.,       & !< flag parameter for unscheduled radiation calls
174                lw_radiation = .TRUE.,                & !< flag parameter indicating whether longwave radiation shall be calculated
175                radiation = .FALSE.,                  & !< flag parameter indicating whether the radiation model is used
176                sun_up    = .TRUE.,                   & !< flag parameter indicating whether the sun is up or down
177                sw_radiation = .TRUE.                   !< flag parameter indicing whether shortwave radiation shall be calculated
[1585]178
[1496]179
[1691]180    REAL(wp), PARAMETER :: d_seconds_hour  = 0.000277777777778_wp,  & !< inverse of seconds per hour (1/3600)
181                           d_hours_day    = 0.0416666666667_wp,     & !< inverse of hours per day (1/24)
182                           sigma_sb       = 5.67037321E-8_wp,       & !< Stefan-Boltzmann constant
183                           solar_constant = 1368.0_wp                 !< solar constant at top of atmosphere
[1585]184
[1691]185    REAL(wp) :: albedo = 9999999.9_wp,           & !< NAMELIST alpha
186                albedo_lw_dif = 9999999.9_wp,    & !< NAMELIST aldif
187                albedo_lw_dir = 9999999.9_wp,    & !< NAMELIST aldir
188                albedo_sw_dif = 9999999.9_wp,    & !< NAMELIST asdif
189                albedo_sw_dir = 9999999.9_wp,    & !< NAMELIST asdir
190                decl_1,                          & !< declination coef. 1
191                decl_2,                          & !< declination coef. 2
192                decl_3,                          & !< declination coef. 3
193                dt_radiation = 0.0_wp,           & !< radiation model timestep
194                emissivity = 0.98_wp,            & !< NAMELIST surface emissivity
195                lambda = 0.0_wp,                 & !< longitude in degrees
196                lon = 0.0_wp,                    & !< longitude in radians
197                lat = 0.0_wp,                    & !< latitude in radians
198                net_radiation = 0.0_wp,          & !< net radiation at surface
199                skip_time_do_radiation = 0.0_wp, & !< Radiation model is not called before this time
200                sky_trans,                       & !< sky transmissivity
201                time_radiation = 0.0_wp,         & !< time since last call of radiation code
202                time_utc,                        & !< current time in UTC
203                time_utc_init = 43200.0_wp         !< UTC time at model start (noon)
204
[1682]205    REAL(wp), DIMENSION(0:0) ::  zenith        !< solar zenith angle
[1585]206
[1496]207    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: &
[1682]208                alpha,                       & !< surface broadband albedo (used for clear-sky scheme)
[1709]209                rad_lw_out_change_0,         & !< change in LW out due to change in surface temperature
[1682]210                rad_net,                     & !< net radiation at the surface
211                rad_net_av                     !< average of rad_net
[1496]212
[1585]213!
214!-- Land surface albedos for solar zenith angle of 60° after Briegleb (1992)     
215!-- (shortwave, longwave, broadband):   sw,      lw,      bb,
[1788]216    REAL(wp), DIMENSION(0:2,1:17), PARAMETER :: albedo_pars = RESHAPE( (/& 
[1585]217                                   0.06_wp, 0.06_wp, 0.06_wp,            & !  1
218                                   0.09_wp, 0.28_wp, 0.19_wp,            & !  2
219                                   0.11_wp, 0.33_wp, 0.23_wp,            & !  3
220                                   0.11_wp, 0.33_wp, 0.23_wp,            & !  4
221                                   0.14_wp, 0.34_wp, 0.25_wp,            & !  5
222                                   0.06_wp, 0.22_wp, 0.14_wp,            & !  6
223                                   0.06_wp, 0.27_wp, 0.17_wp,            & !  7
224                                   0.06_wp, 0.31_wp, 0.19_wp,            & !  8
225                                   0.06_wp, 0.22_wp, 0.14_wp,            & !  9
226                                   0.06_wp, 0.28_wp, 0.18_wp,            & ! 10
227                                   0.35_wp, 0.51_wp, 0.43_wp,            & ! 11
228                                   0.24_wp, 0.40_wp, 0.32_wp,            & ! 12
229                                   0.10_wp, 0.27_wp, 0.19_wp,            & ! 13
230                                   0.90_wp, 0.65_wp, 0.77_wp,            & ! 14
[1587]231                                   0.90_wp, 0.65_wp, 0.77_wp,            & ! 15
[1788]232                                   0.95_wp, 0.70_wp, 0.82_wp,            & ! 16
233                                   0.08_wp, 0.08_wp, 0.08_wp             & ! 17
234                                 /), (/ 3, 17 /) )
[1496]235
[1585]236    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET :: &
[1691]237                        rad_lw_cs_hr,                  & !< longwave clear sky radiation heating rate (K/s)
238                        rad_lw_cs_hr_av,               & !< average of rad_lw_cs_hr
239                        rad_lw_hr,                     & !< longwave radiation heating rate (K/s)
240                        rad_lw_hr_av,                  & !< average of rad_sw_hr
241                        rad_lw_in,                     & !< incoming longwave radiation (W/m2)
242                        rad_lw_in_av,                  & !< average of rad_lw_in
243                        rad_lw_out,                    & !< outgoing longwave radiation (W/m2)
244                        rad_lw_out_av,                 & !< average of rad_lw_out
245                        rad_sw_cs_hr,                  & !< shortwave clear sky radiation heating rate (K/s)
246                        rad_sw_cs_hr_av,               & !< average of rad_sw_cs_hr
247                        rad_sw_hr,                     & !< shortwave radiation heating rate (K/s)
248                        rad_sw_hr_av,                  & !< average of rad_sw_hr
[1682]249                        rad_sw_in,                     & !< incoming shortwave radiation (W/m2)
250                        rad_sw_in_av,                  & !< average of rad_sw_in
251                        rad_sw_out,                    & !< outgoing shortwave radiation (W/m2)
[1691]252                        rad_sw_out_av                    !< average of rad_sw_out
[1585]253
[1691]254
[1585]255!
256!-- Variables and parameters used in RRTMG only
257#if defined ( __rrtmg )
[1682]258    CHARACTER(LEN=12) :: rrtm_input_file = "RAD_SND_DATA" !< name of the NetCDF input file (sounding data)
[1585]259
260
261!
262!-- Flag parameters for RRTMGS (should not be changed)
[1682]263    INTEGER(iwp), PARAMETER :: rrtm_inflglw  = 2, & !< flag for lw cloud optical properties (0,1,2)
264                               rrtm_iceflglw = 0, & !< flag for lw ice particle specifications (0,1,2,3)
265                               rrtm_liqflglw = 1, & !< flag for lw liquid droplet specifications
266                               rrtm_inflgsw  = 2, & !< flag for sw cloud optical properties (0,1,2)
267                               rrtm_iceflgsw = 0, & !< flag for sw ice particle specifications (0,1,2,3)
268                               rrtm_liqflgsw = 1    !< flag for sw liquid droplet specifications
[1585]269
270!
271!-- The following variables should be only changed with care, as this will
272!-- require further setting of some variables, which is currently not
273!-- implemented (aerosols, ice phase).
[1682]274    INTEGER(iwp) :: nzt_rad,           & !< upper vertical limit for radiation calculations
275                    rrtm_icld = 0,     & !< cloud flag (0: clear sky column, 1: cloudy column)
276                    rrtm_iaer = 0,     & !< aerosol option flag (0: no aerosol layers, for lw only: 6 (requires setting of rrtm_sw_ecaer), 10: one or more aerosol layers (not implemented)
[1691]277                    rrtm_idrv = 1        !< longwave upward flux calculation option (0,1)
[1585]278
[1788]279    INTEGER(iwp) :: nc_stat !< local variable for storin the result of netCDF calls for error message handling
280
[1682]281    LOGICAL :: snd_exists = .FALSE.      !< flag parameter to check whether a user-defined input files exists
[1585]282
[1691]283    REAL(wp), PARAMETER :: mol_mass_air_d_wv = 1.607793_wp !< molecular weight dry air / water vapor
[1585]284
[1682]285    REAL(wp), DIMENSION(:), ALLOCATABLE :: hyp_snd,     & !< hypostatic pressure from sounding data (hPa)
286                                           q_snd,       & !< specific humidity from sounding data (kg/kg) - dummy at the moment
287                                           rrtm_tsfc,   & !< dummy array for storing surface temperature
288                                           t_snd          !< actual temperature from sounding data (hPa)
[1585]289
[1691]290    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: aldif,          & !< longwave diffuse albedo solar angle of 60°
291                                             aldir,          & !< longwave direct albedo solar angle of 60°
292                                             asdif,          & !< shortwave diffuse albedo solar angle of 60°
293                                             asdir,          & !< shortwave direct albedo solar angle of 60°
294                                             rrtm_ccl4vmr,   & !< CCL4 volume mixing ratio (g/mol)
295                                             rrtm_cfc11vmr,  & !< CFC11 volume mixing ratio (g/mol)
296                                             rrtm_cfc12vmr,  & !< CFC12 volume mixing ratio (g/mol)
297                                             rrtm_cfc22vmr,  & !< CFC22 volume mixing ratio (g/mol)
298                                             rrtm_ch4vmr,    & !< CH4 volume mixing ratio
299                                             rrtm_cicewp,    & !< in-cloud ice water path (g/m²)
300                                             rrtm_cldfr,     & !< cloud fraction (0,1)
301                                             rrtm_cliqwp,    & !< in-cloud liquid water path (g/m²)
302                                             rrtm_co2vmr,    & !< CO2 volume mixing ratio (g/mol)
303                                             rrtm_emis,      & !< surface emissivity (0-1)   
304                                             rrtm_h2ovmr,    & !< H2O volume mixing ratio
305                                             rrtm_n2ovmr,    & !< N2O volume mixing ratio
306                                             rrtm_o2vmr,     & !< O2 volume mixing ratio
307                                             rrtm_o3vmr,     & !< O3 volume mixing ratio
308                                             rrtm_play,      & !< pressure layers (hPa, zu-grid)
309                                             rrtm_plev,      & !< pressure layers (hPa, zw-grid)
310                                             rrtm_reice,     & !< cloud ice effective radius (microns)
311                                             rrtm_reliq,     & !< cloud water drop effective radius (microns)
312                                             rrtm_tlay,      & !< actual temperature (K, zu-grid)
313                                             rrtm_tlev,      & !< actual temperature (K, zw-grid)
314                                             rrtm_lwdflx,    & !< RRTM output of incoming longwave radiation flux (W/m2)
315                                             rrtm_lwdflxc,   & !< RRTM output of outgoing clear sky longwave radiation flux (W/m2)
316                                             rrtm_lwuflx,    & !< RRTM output of outgoing longwave radiation flux (W/m2)
317                                             rrtm_lwuflxc,   & !< RRTM output of incoming clear sky longwave radiation flux (W/m2)
318                                             rrtm_lwuflx_dt, & !< RRTM output of incoming clear sky longwave radiation flux (W/m2)
319                                             rrtm_lwuflxc_dt,& !< RRTM output of outgoing clear sky longwave radiation flux (W/m2)
320                                             rrtm_lwhr,      & !< RRTM output of longwave radiation heating rate (K/d)
321                                             rrtm_lwhrc,     & !< RRTM output of incoming longwave clear sky radiation heating rate (K/d)
322                                             rrtm_swdflx,    & !< RRTM output of incoming shortwave radiation flux (W/m2)
323                                             rrtm_swdflxc,   & !< RRTM output of outgoing clear sky shortwave radiation flux (W/m2)
324                                             rrtm_swuflx,    & !< RRTM output of outgoing shortwave radiation flux (W/m2)
325                                             rrtm_swuflxc,   & !< RRTM output of incoming clear sky shortwave radiation flux (W/m2)
326                                             rrtm_swhr,      & !< RRTM output of shortwave radiation heating rate (K/d)
327                                             rrtm_swhrc        !< RRTM output of incoming shortwave clear sky radiation heating rate (K/d)
[1585]328
329!
330!-- Definition of arrays that are currently not used for calling RRTMG (due to setting of flag parameters)
[1682]331    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rad_lw_cs_in,   & !< incoming clear sky longwave radiation (W/m2) (not used)
332                                                rad_lw_cs_out,  & !< outgoing clear sky longwave radiation (W/m2) (not used)
333                                                rad_sw_cs_in,   & !< incoming clear sky shortwave radiation (W/m2) (not used)
334                                                rad_sw_cs_out,  & !< outgoing clear sky shortwave radiation (W/m2) (not used)
335                                                rrtm_aldif,     & !< surface albedo for longwave diffuse radiation
336                                                rrtm_aldir,     & !< surface albedo for longwave direct radiation
337                                                rrtm_asdif,     & !< surface albedo for shortwave diffuse radiation
338                                                rrtm_asdir,     & !< surface albedo for shortwave direct radiation
339                                                rrtm_lw_tauaer, & !< lw aerosol optical depth
340                                                rrtm_lw_taucld, & !< lw in-cloud optical depth
341                                                rrtm_sw_taucld, & !< sw in-cloud optical depth
342                                                rrtm_sw_ssacld, & !< sw in-cloud single scattering albedo
343                                                rrtm_sw_asmcld, & !< sw in-cloud asymmetry parameter
344                                                rrtm_sw_fsfcld, & !< sw in-cloud forward scattering fraction
345                                                rrtm_sw_tauaer, & !< sw aerosol optical depth
346                                                rrtm_sw_ssaaer, & !< sw aerosol single scattering albedo
347                                                rrtm_sw_asmaer, & !< sw aerosol asymmetry parameter
348                                                rrtm_sw_ecaer     !< sw aerosol optical detph at 0.55 microns (rrtm_iaer = 6 only)
[1691]349
[1585]350#endif
351
[1826]352    INTERFACE radiation_check_data_output
353       MODULE PROCEDURE radiation_check_data_output
354    END INTERFACE radiation_check_data_output
[1496]355
[1826]356    INTERFACE radiation_check_data_output_pr
357       MODULE PROCEDURE radiation_check_data_output_pr
358    END INTERFACE radiation_check_data_output_pr
359 
360    INTERFACE radiation_check_parameters
361       MODULE PROCEDURE radiation_check_parameters
362    END INTERFACE radiation_check_parameters
363 
[1551]364    INTERFACE radiation_clearsky
365       MODULE PROCEDURE radiation_clearsky
366    END INTERFACE radiation_clearsky
[1826]367 
368    INTERFACE radiation_header
369       MODULE PROCEDURE radiation_header
370    END INTERFACE radiation_header 
371 
372    INTERFACE radiation_init
373       MODULE PROCEDURE radiation_init
374    END INTERFACE radiation_init
[1496]375
[1826]376    INTERFACE radiation_parin
377       MODULE PROCEDURE radiation_parin
378    END INTERFACE radiation_parin
379   
[1585]380    INTERFACE radiation_rrtmg
381       MODULE PROCEDURE radiation_rrtmg
382    END INTERFACE radiation_rrtmg
[1551]383
[1585]384    INTERFACE radiation_tendency
385       MODULE PROCEDURE radiation_tendency
386       MODULE PROCEDURE radiation_tendency_ij
387    END INTERFACE radiation_tendency
[1551]388
[1496]389    SAVE
390
391    PRIVATE
392
[1826]393!
394!-- Public functions
395    PUBLIC radiation_check_data_output, radiation_check_data_output_pr,        &
396           radiation_check_parameters, radiation_clearsky, radiation_header,   &
397           radiation_init, radiation_parin, radiation_rrtmg, radiation_tendency
398   
399!
400!-- Public variables and constants
401    PUBLIC dots_rad, dt_radiation, force_radiation_call,                       &
402           rad_net, rad_net_av, radiation, radiation_scheme, rad_lw_in,        &
403           rad_lw_in_av, rad_lw_out, rad_lw_out_av, rad_lw_out_change_0,       &
404           rad_lw_cs_hr, rad_lw_cs_hr_av, rad_lw_hr, rad_lw_hr_av, rad_sw_in,  &
405           rad_sw_in_av, rad_sw_out, rad_sw_out_av, rad_sw_cs_hr,              &
406           rad_sw_cs_hr_av, rad_sw_hr, rad_sw_hr_av, sigma_sb,                 &
407           skip_time_do_radiation, time_radiation, unscheduled_radiation_calls
[1496]408
[1691]409
[1585]410#if defined ( __rrtmg )
[1709]411    PUBLIC rrtm_aldif, rrtm_aldir, rrtm_asdif, rrtm_asdir, rrtm_idrv
[1585]412#endif
[1496]413
414 CONTAINS
415
416!------------------------------------------------------------------------------!
417! Description:
418! ------------
[1826]419!> Check data output for radiation model
420!------------------------------------------------------------------------------!
421    SUBROUTINE radiation_check_data_output( var, unit, i, ilen, k )
422 
423 
424       USE control_parameters,                                                 &
425           ONLY: data_output, message_string
426
427       IMPLICIT NONE
428
429       CHARACTER (LEN=*) ::  unit     !<
430       CHARACTER (LEN=*) ::  var      !<
431
432       INTEGER(iwp) :: i
433       INTEGER(iwp) :: ilen
434       INTEGER(iwp) :: k
435
436       SELECT CASE ( TRIM( var ) )
437
438         CASE ( 'rad_lw_in', 'rad_lw_out', 'rad_lw_cs_hr', 'rad_lw_hr',       &
439                 'rad_sw_in', 'rad_sw_out', 'rad_sw_cs_hr', 'rad_sw_hr' )
440             IF (  .NOT.  radiation  .OR.  radiation_scheme /= 'rrtmg' )  THEN
441                message_string = '"output of "' // TRIM( var ) // '" requi' // &
442                                 'res radiation = .TRUE. and ' //              &
443                                 'radiation_scheme = "rrtmg"'
444                CALL message( 'check_parameters', 'PA0406', 1, 2, 0, 6, 0 )
445             ENDIF
446             unit = 'W/m2'     
447
448          CASE ( 'rad_net*', 'rrtm_aldif*', 'rrtm_aldir*', 'rrtm_asdif*',      &
449                 'rrtm_asdir*' )
450             IF ( k == 0  .OR.  data_output(i)(ilen-2:ilen) /= '_xy' )  THEN
451                message_string = 'illegal value for data_output: "' //         &
452                                 TRIM( var ) // '" & only 2d-horizontal ' //   &
453                                 'cross sections are allowed for this value'
454                CALL message( 'check_parameters', 'PA0111', 1, 2, 0, 6, 0 )
455             ENDIF
456             IF (  .NOT.  radiation  .OR.  radiation_scheme /= "rrtmg" )  THEN
457                IF ( TRIM( var ) == 'rrtm_aldif*'  .OR.                        &
458                     TRIM( var ) == 'rrtm_aldir*'  .OR.                        &
459                     TRIM( var ) == 'rrtm_asdif*'  .OR.                        &
460                     TRIM( var ) == 'rrtm_asdir*'      )                       &
461                THEN
462                   message_string = 'output of "' // TRIM( var ) // '" require'&
463                                    // 's radiation = .TRUE. and radiation_sch'&
464                                    // 'eme = "rrtmg"'
465                   CALL message( 'check_parameters', 'PA0409', 1, 2, 0, 6, 0 )
466                ENDIF
467             ENDIF
468
469             IF ( TRIM( var ) == 'rad_net*'      ) unit = 'W/m2'   
470             IF ( TRIM( var ) == 'rrtm_aldif*'   ) unit = ''   
471             IF ( TRIM( var ) == 'rrtm_aldir*'   ) unit = '' 
472             IF ( TRIM( var ) == 'rrtm_asdif*'   ) unit = '' 
473             IF ( TRIM( var ) == 'rrtm_asdir*'   ) unit = '' 
474
475          CASE DEFAULT
476             unit = 'illegal'
477
478       END SELECT
479
480
481    END SUBROUTINE radiation_check_data_output
482
483!------------------------------------------------------------------------------!
484! Description:
485! ------------
486!> Check data output of profiles for radiation model
487!------------------------------------------------------------------------------! 
488    SUBROUTINE radiation_check_data_output_pr( variable, var_count, unit, dopr_unit )
489 
490       USE arrays_3d,                                                          &
491           ONLY: zu
492
493       USE control_parameters,                                                 &
494           ONLY: data_output_pr, message_string
495
496       USE indices
497
498       USE profil_parameter
499
500       USE statistics
501
502       IMPLICIT NONE
503   
504       CHARACTER (LEN=*) ::  unit      !<
505       CHARACTER (LEN=*) ::  variable  !<
506       CHARACTER (LEN=*) ::  dopr_unit !< local value of dopr_unit
507 
508       INTEGER(iwp) ::  user_pr_index !<
509       INTEGER(iwp) ::  var_count     !<
510
511       SELECT CASE ( TRIM( variable ) )
512       
513         CASE ( 'rad_net' )
514             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme == 'constant' )&
515             THEN
516                message_string = 'data_output_pr = ' //                        &
517                                 TRIM( data_output_pr(var_count) ) // ' is' // &
518                                 'not available for radiation = .FALSE. or ' //&
519                                 'radiation_scheme = "constant"'
520                CALL message( 'check_parameters', 'PA0408', 1, 2, 0, 6, 0 )
521             ELSE
522                dopr_index(var_count) = 101
523                dopr_unit  = 'W/m2'
524                hom(:,2,101,:)  = SPREAD( zw, 2, statistic_regions+1 )
525                unit = dopr_unit
526             ENDIF
527
528          CASE ( 'rad_lw_in' )
529             IF ( (  .NOT.  radiation)  .OR.  radiation_scheme == 'constant' ) &
530             THEN
531                message_string = 'data_output_pr = ' //                        &
532                                 TRIM( data_output_pr(var_count) ) // ' is' // &
533                                 'not available for radiation = .FALSE. or ' //&
534                                 'radiation_scheme = "constant"'
535                CALL message( 'check_parameters', 'PA0408', 1, 2, 0, 6, 0 )
536             ELSE
537                dopr_index(var_count) = 102
538                dopr_unit  = 'W/m2'
539                hom(:,2,102,:)  = SPREAD( zw, 2, statistic_regions+1 )
540                unit = dopr_unit 
541             ENDIF
542
543          CASE ( 'rad_lw_out' )
544             IF ( (  .NOT. radiation )  .OR.  radiation_scheme == 'constant' ) &
545             THEN
546                message_string = 'data_output_pr = ' //                        &
547                                 TRIM( data_output_pr(var_count) ) // ' is' // &
548                                 'not available for radiation = .FALSE. or ' //&
549                                 'radiation_scheme = "constant"'
550                CALL message( 'check_parameters', 'PA0408', 1, 2, 0, 6, 0 )
551             ELSE
552                dopr_index(var_count) = 103
553                dopr_unit  = 'W/m2'
554                hom(:,2,103,:)  = SPREAD( zw, 2, statistic_regions+1 )
555                unit = dopr_unit   
556             ENDIF
557
558          CASE ( 'rad_sw_in' )
559             IF ( (  .NOT. radiation )  .OR.  radiation_scheme == 'constant' ) &
560             THEN
561                message_string = 'data_output_pr = ' //                        &
562                                 TRIM( data_output_pr(var_count) ) // ' is' // &
563                                 'not available for radiation = .FALSE. or ' //&
564                                 'radiation_scheme = "constant"'
565                CALL message( 'check_parameters', 'PA0408', 1, 2, 0, 6, 0 )
566             ELSE
567                dopr_index(var_count) = 104
568                dopr_unit  = 'W/m2'
569                hom(:,2,104,:)  = SPREAD( zw, 2, statistic_regions+1 )
570                unit = dopr_unit
571             ENDIF
572
573          CASE ( 'rad_sw_out')
574             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme == 'constant' )&
575             THEN
576                message_string = 'data_output_pr = ' //                        &
577                                 TRIM( data_output_pr(var_count) ) // ' is' // &
578                                 'not available for radiation = .FALSE. or ' //&
579                                 'radiation_scheme = "constant"'
580                CALL message( 'check_parameters', 'PA0408', 1, 2, 0, 6, 0 )
581             ELSE
582                dopr_index(var_count) = 105
583                dopr_unit  = 'W/m2'
584                hom(:,2,105,:)  = SPREAD( zw, 2, statistic_regions+1 )
585                unit = dopr_unit
586             ENDIF
587
588          CASE ( 'rad_lw_cs_hr' )
589             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme /= 'rrtmg' )   &
590             THEN
591                message_string = 'data_output_pr = ' //                        &
592                                 TRIM( data_output_pr(var_count) ) // ' is' // &
593                                 'not available for radiation = .FALSE. or ' //&
594                                 'radiation_scheme /= "rrtmg"'
595                CALL message( 'check_parameters', 'PA0413', 1, 2, 0, 6, 0 )
596             ELSE
597                dopr_index(var_count) = 106
598                dopr_unit  = 'K/h'
599                hom(:,2,106,:)  = SPREAD( zu, 2, statistic_regions+1 )
600                unit = dopr_unit
601             ENDIF
602
603          CASE ( 'rad_lw_hr' )
604             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme /= 'rrtmg' )   &
605             THEN
606                message_string = 'data_output_pr = ' //                        &
607                                 TRIM( data_output_pr(var_count) ) // ' is' // &
608                                 'not available for radiation = .FALSE. or ' //&
609                                 'radiation_scheme /= "rrtmg"'
610                CALL message( 'check_parameters', 'PA0413', 1, 2, 0, 6, 0 )
611             ELSE
612                dopr_index(var_count) = 107
613                dopr_unit  = 'K/h'
614                hom(:,2,107,:)  = SPREAD( zu, 2, statistic_regions+1 )
615                unit = dopr_unit
616             ENDIF
617
618          CASE ( 'rad_sw_cs_hr' )
619             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme /= 'rrtmg' )   &
620             THEN
621                message_string = 'data_output_pr = ' //                        &
622                                 TRIM( data_output_pr(var_count) ) // ' is' // &
623                                 'not available for radiation = .FALSE. or ' //&
624                                 'radiation_scheme /= "rrtmg"'
625                CALL message( 'check_parameters', 'PA0413', 1, 2, 0, 6, 0 )
626             ELSE
627                dopr_index(var_count) = 108
628                dopr_unit  = 'K/h'
629                hom(:,2,108,:)  = SPREAD( zu, 2, statistic_regions+1 )
630                unit = dopr_unit
631             ENDIF
632
633          CASE ( 'rad_sw_hr' )
634             IF ( (  .NOT.  radiation )  .OR.  radiation_scheme /= 'rrtmg' )   &
635             THEN
636                message_string = 'data_output_pr = ' //                        &
637                                 TRIM( data_output_pr(var_count) ) // ' is' // &
638                                 'not available for radiation = .FALSE. or ' //&
639                                 'radiation_scheme /= "rrtmg"'
640                CALL message( 'check_parameters', 'PA0413', 1, 2, 0, 6, 0 )
641             ELSE
642                dopr_index(var_count) = 109
643                dopr_unit  = 'K/h'
644                hom(:,2,109,:)  = SPREAD( zu, 2, statistic_regions+1 )
645                unit = dopr_unit
646             ENDIF
647
648
649          CASE DEFAULT
650             unit = 'illegal'
651
652       END SELECT
653
654
655    END SUBROUTINE radiation_check_data_output_pr
656 
657 
658!------------------------------------------------------------------------------!
659! Description:
660! ------------
661!> Check parameters routine for radiation model
662!------------------------------------------------------------------------------!
663    SUBROUTINE radiation_check_parameters
664
665       USE control_parameters,                                                 &
666           ONLY: message_string, topography
667                 
668   
669       IMPLICIT NONE
670
671       IF ( radiation_scheme /= 'constant'   .AND.                             &
672            radiation_scheme /= 'clear-sky'  .AND.                             &
673            radiation_scheme /= 'rrtmg' )  THEN
674          message_string = 'unknown radiation_scheme = '//                     &
675                           TRIM( radiation_scheme )
676          CALL message( 'check_parameters', 'PA0405', 1, 2, 0, 6, 0 )
677       ELSEIF ( radiation_scheme == 'rrtmg' )  THEN
678#if ! defined ( __rrtmg )
679          message_string = 'radiation_scheme = "rrtmg" requires ' //           & 
680                           'compilation of PALM with pre-processor ' //        &
681                           'directive -D__rrtmg'
682          CALL message( 'check_parameters', 'PA0407', 1, 2, 0, 6, 0 )
683#endif
684#if defined ( __rrtmg ) && ! defined( __netcdf )
685          message_string = 'radiation_scheme = "rrtmg" requires ' //           & 
686                           'the use of NetCDF (preprocessor directive ' //     &
687                           '-D__netcdf'
688          CALL message( 'check_parameters', 'PA0412', 1, 2, 0, 6, 0 )
689#endif
690
691       ENDIF
692
693       IF ( albedo_type == 0  .AND.  albedo == 9999999.9_wp  .AND.             &
694            radiation_scheme == 'clear-sky')  THEN
695          message_string = 'radiation_scheme = "clear-sky" in combination' //  & 
696                           'with albedo_type = 0 requires setting of albedo'// &
697                           ' /= 9999999.9'
698          CALL message( 'check_parameters', 'PA0410', 1, 2, 0, 6, 0 )
699       ENDIF
700
701       IF ( albedo_type == 0  .AND.  radiation_scheme == 'rrtmg'  .AND.        &
702          (    albedo_lw_dif == 9999999.9_wp .OR. albedo_lw_dir == 9999999.9_wp&
703          .OR. albedo_sw_dif == 9999999.9_wp .OR. albedo_sw_dir == 9999999.9_wp& 
704          ) ) THEN
705          message_string = 'radiation_scheme = "rrtmg" in combination' //      & 
706                           'with albedo_type = 0 requires setting of ' //      &
707                           'albedo_lw_dif /= 9999999.9' //                     &
708                           'albedo_lw_dir /= 9999999.9' //                     &
709                           'albedo_sw_dif /= 9999999.9 and' //                 &
710                           'albedo_sw_dir /= 9999999.9'
711          CALL message( 'check_parameters', 'PA0411', 1, 2, 0, 6, 0 )
712       ENDIF
713
714       IF ( topography /= 'flat' )  THEN
715          message_string = 'radiation scheme cannot be used ' //               & 
716                           'in combination with  topography /= "flat"'
717          CALL message( 'check_parameters', 'PA0414', 1, 2, 0, 6, 0 )
718       ENDIF
719 
720    END SUBROUTINE radiation_check_parameters 
721 
722 
723!------------------------------------------------------------------------------!
724! Description:
725! ------------
[1682]726!> Initialization of the radiation model
[1496]727!------------------------------------------------------------------------------!
[1826]728    SUBROUTINE radiation_init
[1496]729   
730       IMPLICIT NONE
731
[1585]732!
733!--    Allocate array for storing the surface net radiation
734       IF ( .NOT. ALLOCATED ( rad_net ) )  THEN
735          ALLOCATE ( rad_net(nysg:nyng,nxlg:nxrg) )
736          rad_net = 0.0_wp
737       ENDIF
[1496]738
739!
[1709]740!--    Allocate array for storing the surface net radiation
741       IF ( .NOT. ALLOCATED ( rad_lw_out_change_0 ) )  THEN
742          ALLOCATE ( rad_lw_out_change_0(nysg:nyng,nxlg:nxrg) )
743          rad_lw_out_change_0 = 0.0_wp
744       ENDIF
745
746!
[1551]747!--    Fix net radiation in case of radiation_scheme = 'constant'
[1585]748       IF ( radiation_scheme == 'constant' )  THEN
[1551]749          rad_net = net_radiation
[1585]750          radiation = .FALSE.
[1551]751!
[1585]752!--    Calculate orbital constants
753       ELSE
[1551]754          decl_1 = SIN(23.45_wp * pi / 180.0_wp)
755          decl_2 = 2.0_wp * pi / 365.0_wp
756          decl_3 = decl_2 * 81.0_wp
[1585]757          lat    = phi * pi / 180.0_wp
758          lon    = lambda * pi / 180.0_wp
759       ENDIF
760
761
762       IF ( radiation_scheme == 'clear-sky' )  THEN
763
764          ALLOCATE ( alpha(nysg:nyng,nxlg:nxrg) )
765
766          IF ( .NOT. ALLOCATED ( rad_sw_in ) )  THEN
767             ALLOCATE ( rad_sw_in(0:0,nysg:nyng,nxlg:nxrg) )
768          ENDIF
769          IF ( .NOT. ALLOCATED ( rad_sw_out ) )  THEN
770             ALLOCATE ( rad_sw_out(0:0,nysg:nyng,nxlg:nxrg) )
771          ENDIF
772
773          IF ( .NOT. ALLOCATED ( rad_sw_in_av ) )  THEN
774             ALLOCATE ( rad_sw_in_av(0:0,nysg:nyng,nxlg:nxrg) )
775          ENDIF
776          IF ( .NOT. ALLOCATED ( rad_sw_out_av ) )  THEN
777             ALLOCATE ( rad_sw_out_av(0:0,nysg:nyng,nxlg:nxrg) )
778          ENDIF
779
780          IF ( .NOT. ALLOCATED ( rad_lw_in ) )  THEN
781             ALLOCATE ( rad_lw_in(0:0,nysg:nyng,nxlg:nxrg) )
782          ENDIF
783          IF ( .NOT. ALLOCATED ( rad_lw_out ) )  THEN
784             ALLOCATE ( rad_lw_out(0:0,nysg:nyng,nxlg:nxrg) )
785          ENDIF
786
787          IF ( .NOT. ALLOCATED ( rad_lw_in_av ) )  THEN
788             ALLOCATE ( rad_lw_in_av(0:0,nysg:nyng,nxlg:nxrg) )
789          ENDIF
790          IF ( .NOT. ALLOCATED ( rad_lw_out_av ) )  THEN
791             ALLOCATE ( rad_lw_out_av(0:0,nysg:nyng,nxlg:nxrg) )
792          ENDIF
793
794          rad_sw_in  = 0.0_wp
795          rad_sw_out = 0.0_wp
796          rad_lw_in  = 0.0_wp
797          rad_lw_out = 0.0_wp
798
[1496]799!
[1585]800!--       Overwrite albedo if manually set in parameter file
801          IF ( albedo_type /= 0 .AND. albedo == 9999999.9_wp )  THEN
802             albedo = albedo_pars(2,albedo_type)
803          ENDIF
804   
805          alpha = albedo
806 
807!
808!--    Initialization actions for RRTMG
809       ELSEIF ( radiation_scheme == 'rrtmg' )  THEN
810#if defined ( __rrtmg )
811!
812!--       Allocate albedos
813          ALLOCATE ( rrtm_aldif(0:0,nysg:nyng,nxlg:nxrg) )
814          ALLOCATE ( rrtm_aldir(0:0,nysg:nyng,nxlg:nxrg) )
815          ALLOCATE ( rrtm_asdif(0:0,nysg:nyng,nxlg:nxrg) )
816          ALLOCATE ( rrtm_asdir(0:0,nysg:nyng,nxlg:nxrg) )
817          ALLOCATE ( aldif(nysg:nyng,nxlg:nxrg) )
818          ALLOCATE ( aldir(nysg:nyng,nxlg:nxrg) )
819          ALLOCATE ( asdif(nysg:nyng,nxlg:nxrg) )
820          ALLOCATE ( asdir(nysg:nyng,nxlg:nxrg) )
821
822          IF ( albedo_type /= 0 )  THEN
823             IF ( albedo_lw_dif == 9999999.9_wp )  THEN
824                albedo_lw_dif = albedo_pars(0,albedo_type)
825                albedo_lw_dir = albedo_lw_dif
826             ENDIF
827             IF ( albedo_sw_dif == 9999999.9_wp )  THEN
828                albedo_sw_dif = albedo_pars(1,albedo_type)
829                albedo_sw_dir = albedo_sw_dif
830             ENDIF
831          ENDIF
832
833          aldif(:,:) = albedo_lw_dif
834          aldir(:,:) = albedo_lw_dir
835          asdif(:,:) = albedo_sw_dif
836          asdir(:,:) = albedo_sw_dir
837!
838!--       Calculate initial values of current (cosine of) the zenith angle and
839!--       whether the sun is up
840          CALL calc_zenith     
841!
842!--       Calculate initial surface albedo
843          IF ( .NOT. constant_albedo )  THEN
844             CALL calc_albedo
845          ELSE
846             rrtm_aldif(0,:,:) = aldif(:,:)
847             rrtm_aldir(0,:,:) = aldir(:,:)
848             rrtm_asdif(0,:,:) = asdif(:,:) 
849             rrtm_asdir(0,:,:) = asdir(:,:)   
850          ENDIF
851
852!
853!--       Allocate surface emissivity
854          ALLOCATE ( rrtm_emis(0:0,1:nbndlw+1) )
855          rrtm_emis = emissivity
856
857!
858!--       Allocate 3d arrays of radiative fluxes and heating rates
859          IF ( .NOT. ALLOCATED ( rad_sw_in ) )  THEN
860             ALLOCATE ( rad_sw_in(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
861             rad_sw_in = 0.0_wp
862          ENDIF
863
864          IF ( .NOT. ALLOCATED ( rad_sw_in_av ) )  THEN
865             ALLOCATE ( rad_sw_in_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
866          ENDIF
867
868          IF ( .NOT. ALLOCATED ( rad_sw_out ) )  THEN
869             ALLOCATE ( rad_sw_out(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1691]870             rad_sw_out = 0.0_wp
[1585]871          ENDIF
872
873          IF ( .NOT. ALLOCATED ( rad_sw_out_av ) )  THEN
874             ALLOCATE ( rad_sw_out_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
875          ENDIF
876
[1691]877          IF ( .NOT. ALLOCATED ( rad_sw_hr ) )  THEN
878             ALLOCATE ( rad_sw_hr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
879             rad_sw_hr = 0.0_wp
880          ENDIF
[1585]881
[1691]882          IF ( .NOT. ALLOCATED ( rad_sw_hr_av ) )  THEN
883             ALLOCATE ( rad_sw_hr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
884             rad_sw_hr_av = 0.0_wp
885          ENDIF
886
887          IF ( .NOT. ALLOCATED ( rad_sw_cs_hr ) )  THEN
888             ALLOCATE ( rad_sw_cs_hr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
889             rad_sw_cs_hr = 0.0_wp
890          ENDIF
891
892          IF ( .NOT. ALLOCATED ( rad_sw_cs_hr_av ) )  THEN
893             ALLOCATE ( rad_sw_cs_hr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
894             rad_sw_cs_hr_av = 0.0_wp
895          ENDIF
896
[1585]897          IF ( .NOT. ALLOCATED ( rad_lw_in ) )  THEN
898             ALLOCATE ( rad_lw_in(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
899             rad_lw_in     = 0.0_wp
900          ENDIF
901
902          IF ( .NOT. ALLOCATED ( rad_lw_in_av ) )  THEN
903             ALLOCATE ( rad_lw_in_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
904          ENDIF
905
906          IF ( .NOT. ALLOCATED ( rad_lw_out ) )  THEN
907             ALLOCATE ( rad_lw_out(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
908            rad_lw_out    = 0.0_wp
909          ENDIF
910
911          IF ( .NOT. ALLOCATED ( rad_lw_out_av ) )  THEN
912             ALLOCATE ( rad_lw_out_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
913          ENDIF
914
[1691]915          IF ( .NOT. ALLOCATED ( rad_lw_hr ) )  THEN
916             ALLOCATE ( rad_lw_hr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
917             rad_lw_hr = 0.0_wp
918          ENDIF
919
920          IF ( .NOT. ALLOCATED ( rad_lw_hr_av ) )  THEN
921             ALLOCATE ( rad_lw_hr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
922             rad_lw_hr_av = 0.0_wp
923          ENDIF
924
925          IF ( .NOT. ALLOCATED ( rad_lw_cs_hr ) )  THEN
926             ALLOCATE ( rad_lw_cs_hr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
927             rad_lw_cs_hr = 0.0_wp
928          ENDIF
929
930          IF ( .NOT. ALLOCATED ( rad_lw_cs_hr_av ) )  THEN
931             ALLOCATE ( rad_lw_cs_hr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
932             rad_lw_cs_hr_av = 0.0_wp
933          ENDIF
934
935          ALLOCATE ( rad_sw_cs_in(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
936          ALLOCATE ( rad_sw_cs_out(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1585]937          rad_sw_cs_in  = 0.0_wp
938          rad_sw_cs_out = 0.0_wp
939
[1691]940          ALLOCATE ( rad_lw_cs_in(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
941          ALLOCATE ( rad_lw_cs_out(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1585]942          rad_lw_cs_in  = 0.0_wp
943          rad_lw_cs_out = 0.0_wp
944
945!
946!--       Allocate dummy array for storing surface temperature
947          ALLOCATE ( rrtm_tsfc(1) )
948
949!
950!--       Initialize RRTMG
951          IF ( lw_radiation )  CALL rrtmg_lw_ini ( cp )
952          IF ( sw_radiation )  CALL rrtmg_sw_ini ( cp )
953
954!
955!--       Set input files for RRTMG
956          INQUIRE(FILE="RAD_SND_DATA", EXIST=snd_exists) 
957          IF ( .NOT. snd_exists )  THEN
958             rrtm_input_file = "rrtmg_lw.nc"
959          ENDIF
960
961!
962!--       Read vertical layers for RRTMG from sounding data
963!--       The routine provides nzt_rad, hyp_snd(1:nzt_rad),
964!--       t_snd(nzt+2:nzt_rad), rrtm_play(1:nzt_rad), rrtm_plev(1_nzt_rad+1),
965!--       rrtm_tlay(nzt+2:nzt_rad), rrtm_tlev(nzt+2:nzt_rad+1)
966          CALL read_sounding_data
967
968!
969!--       Read trace gas profiles from file. This routine provides
970!--       the rrtm_ arrays (1:nzt_rad+1)
971          CALL read_trace_gas_data
972#endif
[1551]973       ENDIF
[1585]974
[1551]975!
[1585]976!--    Perform user actions if required
977       CALL user_init_radiation
978
979!
980!--    Calculate radiative fluxes at model start
981       IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
982          IF ( radiation_scheme == 'clear-sky' )  THEN
[1709]983             CALL radiation_clearsky
[1585]984          ELSEIF ( radiation_scheme == 'rrtmg' )  THEN
985             CALL radiation_rrtmg
986          ENDIF
987       ENDIF
988
[1496]989       RETURN
990
[1826]991    END SUBROUTINE radiation_init
[1496]992
993
994!------------------------------------------------------------------------------!
995! Description:
996! ------------
[1682]997!> A simple clear sky radiation model
[1496]998!------------------------------------------------------------------------------!
[1551]999    SUBROUTINE radiation_clearsky
[1496]1000
[1585]1001       USE indices,                                                            &
1002           ONLY:  nbgp
1003
[1496]1004       IMPLICIT NONE
1005
[1691]1006       INTEGER(iwp) :: i, j, k   !< loop indices
1007       REAL(wp)     :: exn,   &  !< Exner functions at surface
[1709]1008                       exn1,  &  !< Exner functions at first grid level
1009                       pt1       !< potential temperature at first grid level
[1585]1010
[1496]1011!
[1585]1012!--    Calculate current zenith angle
1013       CALL calc_zenith
1014
1015!
1016!--    Calculate sky transmissivity
1017       sky_trans = 0.6_wp + 0.2_wp * zenith(0)
1018
1019!
1020!--    Calculate value of the Exner function
1021       exn = (surface_pressure / 1000.0_wp )**0.286_wp
1022!
1023!--    Calculate radiation fluxes and net radiation (rad_net) for each grid
1024!--    point
[1709]1025       DO i = nxlg, nxrg
1026          DO j = nysg, nyng
[1585]1027             k = nzb_s_inner(j,i)
[1691]1028
[1709]1029             exn1 = (hyp(k+1) / 100000.0_wp )**0.286_wp
[1691]1030
[1585]1031             rad_sw_in(0,j,i)  = solar_constant * sky_trans * zenith(0)
1032             rad_sw_out(0,j,i) = alpha(j,i) * rad_sw_in(0,j,i)
[1691]1033             rad_lw_out(0,j,i) = emissivity * sigma_sb * (pt(k,j,i) * exn)**4
[1585]1034
[1691]1035             IF ( cloud_physics )  THEN
[1709]1036                pt1 = pt(k+1,j,i) + l_d_cp / exn1 * ql(k+1,j,i)
1037                rad_lw_in(0,j,i)  = 0.8_wp * sigma_sb * (pt1 * exn1)**4
[1691]1038             ELSE
[1709]1039                rad_lw_in(0,j,i)  = 0.8_wp * sigma_sb * (pt(k+1,j,i) * exn1)**4
[1691]1040             ENDIF
1041
1042             rad_net(j,i) = rad_sw_in(0,j,i) - rad_sw_out(0,j,i)               &
1043                            + rad_lw_in(0,j,i) - rad_lw_out(0,j,i)
1044
[1585]1045          ENDDO
1046       ENDDO
1047
1048    END SUBROUTINE radiation_clearsky
1049
1050
1051!------------------------------------------------------------------------------!
1052! Description:
1053! ------------
[1826]1054!> Header output for radiation model
1055!------------------------------------------------------------------------------!
1056    SUBROUTINE radiation_header ( io )
1057
1058
1059       IMPLICIT NONE
1060 
1061       INTEGER(iwp), INTENT(IN) ::  io            !< Unit of the output file
1062   
1063
1064       
1065!
1066!--    Write radiation model header
1067       WRITE( io, 3 )
1068
1069       IF ( radiation_scheme == "constant" )  THEN
1070          WRITE( io, 4 ) net_radiation
1071       ELSEIF ( radiation_scheme == "clear-sky" )  THEN
1072          WRITE( io, 5 )
1073       ELSEIF ( radiation_scheme == "rrtmg" )  THEN
1074          WRITE( io, 6 )
1075          IF ( .NOT. lw_radiation )  WRITE( io, 10 )
1076          IF ( .NOT. sw_radiation )  WRITE( io, 11 )
1077       ENDIF
1078
1079       IF ( albedo_type == 0 )  THEN
1080          WRITE( io, 7 ) albedo
1081       ELSE
1082          WRITE( io, 8 ) TRIM( albedo_type_name(albedo_type) )
1083       ENDIF
1084       IF ( constant_albedo )  THEN
1085          WRITE( io, 9 )
1086       ENDIF
1087       
1088       IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1089          WRITE ( io, 1 )  lambda
1090          WRITE ( io, 2 )  day_init, time_utc_init
1091       ENDIF
1092
1093       WRITE( io, 12 ) dt_radiation
1094 
1095
1096 1 FORMAT ('    Geograph. longitude            :   lambda = ',F4.1,' degr')
1097 2 FORMAT ('    Day of the year at model start :   day_init = ',I3   &
1098            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
1099 3 FORMAT (//' Radiation model information:'/                                  &
1100              ' ----------------------------'/)
1101 4 FORMAT ('    --> Using constant net radiation: net_radiation = ', F6.2,        &
1102           // 'W/m**2')
1103 5 FORMAT ('    --> Simple radiation scheme for clear sky is used (no clouds,',   &
1104                   ' default)')
1105 6 FORMAT ('    --> RRTMG scheme is used')
1106 7 FORMAT (/'    User-specific surface albedo: albedo =', F6.3)
1107 8 FORMAT (/'    Albedo is set for land surface type: ', A)
1108 9 FORMAT (/'    --> Albedo is fixed during the run')
110910 FORMAT (/'    --> Longwave radiation is disabled')
111011 FORMAT (/'    --> Shortwave radiation is disabled.')
111112 FORMAT  ('    Timestep: dt_radiation = ', F6.2, '  s')
1112
1113
1114    END SUBROUTINE radiation_header
1115   
1116
1117!------------------------------------------------------------------------------!
1118! Description:
1119! ------------
1120!> Parin for &radiation_par for radiation model
1121!------------------------------------------------------------------------------!
1122    SUBROUTINE radiation_parin
1123
1124
1125       IMPLICIT NONE
1126
1127       CHARACTER (LEN=80) ::  line  !< dummy string that contains the current line of the parameter file
1128       
1129       NAMELIST /radiation_par/   albedo, albedo_type, albedo_lw_dir,          &
1130                                  albedo_lw_dif, albedo_sw_dir, albedo_sw_dif, &
1131                                  constant_albedo, day_init, dt_radiation,     &
1132                                  lambda, lw_radiation, net_radiation,         &
1133                                  radiation_scheme, skip_time_do_radiation,    &
1134                                  sw_radiation, time_utc_init,                 &
1135                                  unscheduled_radiation_calls
1136       
1137       line = ' '
1138       
1139!
1140!--    Try to find radiation model package
1141       REWIND ( 11 )
1142       line = ' '
1143       DO   WHILE ( INDEX( line, '&radiation_par' ) == 0 )
1144          READ ( 11, '(A)', END=10 )  line
1145       ENDDO
1146       BACKSPACE ( 11 )
1147
1148!
1149!--    Read user-defined namelist
1150       READ ( 11, radiation_par )
1151
1152!
1153!--    Set flag that indicates that the radiation model is switched on
1154       radiation = .TRUE.
1155
1156 10    CONTINUE
1157       
1158
1159    END SUBROUTINE radiation_parin
1160
1161
1162!------------------------------------------------------------------------------!
1163! Description:
1164! ------------
[1682]1165!> Implementation of the RRTMG radiation_scheme
[1585]1166!------------------------------------------------------------------------------!
1167    SUBROUTINE radiation_rrtmg
1168
1169       USE indices,                                                            &
1170           ONLY:  nbgp
1171
1172       USE particle_attributes,                                                &
1173           ONLY:  grid_particles, number_of_particles, particles,              &
1174                  particle_advection_start, prt_count
1175
1176       IMPLICIT NONE
1177
1178#if defined ( __rrtmg )
1179
[1691]1180       INTEGER(iwp) :: i, j, k, n !< loop indices
[1585]1181
[1691]1182       REAL(wp)     ::  s_r2, &   !< weighted sum over all droplets with r^2
1183                        s_r3      !< weighted sum over all droplets with r^3
[1585]1184
1185!
1186!--    Calculate current (cosine of) zenith angle and whether the sun is up
1187       CALL calc_zenith     
1188!
1189!--    Calculate surface albedo
1190       IF ( .NOT. constant_albedo )  THEN
1191          CALL calc_albedo
1192       ENDIF
1193
1194!
1195!--    Prepare input data for RRTMG
1196
1197!
1198!--    In case of large scale forcing with surface data, calculate new pressure
1199!--    profile. nzt_rad might be modified by these calls and all required arrays
1200!--    will then be re-allocated
[1691]1201       IF ( large_scale_forcing  .AND.  lsf_surf )  THEN
[1585]1202          CALL read_sounding_data
1203          CALL read_trace_gas_data
1204       ENDIF
1205!
1206!--    Loop over all grid points
1207       DO i = nxl, nxr
1208          DO j = nys, nyn
1209
1210!
1211!--          Prepare profiles of temperature and H2O volume mixing ratio
[1691]1212             rrtm_tlev(0,nzb+1) = pt(nzb,j,i) * ( surface_pressure             &
1213                                                  / 1000.0_wp )**0.286_wp
[1585]1214
1215             DO k = nzb+1, nzt+1
1216                rrtm_tlay(0,k) = pt(k,j,i) * ( (hyp(k) ) / 100000.0_wp         &
[1691]1217                                 )**0.286_wp + l_d_cp * ql(k,j,i)
1218                rrtm_h2ovmr(0,k) = mol_mass_air_d_wv * (q(k,j,i) - ql(k,j,i))
[1585]1219
1220             ENDDO
1221
1222!
1223!--          Avoid temperature/humidity jumps at the top of the LES domain by
1224!--          linear interpolation from nzt+2 to nzt+7
1225             DO k = nzt+2, nzt+7
1226                rrtm_tlay(0,k) = rrtm_tlay(0,nzt+1)                            &
1227                              + ( rrtm_tlay(0,nzt+8) - rrtm_tlay(0,nzt+1) )    &
1228                              / ( rrtm_play(0,nzt+8) - rrtm_play(0,nzt+1) )    &
1229                              * ( rrtm_play(0,k) - rrtm_play(0,nzt+1) )
1230
1231                rrtm_h2ovmr(0,k) = rrtm_h2ovmr(0,nzt+1)                        &
1232                              + ( rrtm_h2ovmr(0,nzt+8) - rrtm_h2ovmr(0,nzt+1) )&
1233                              / ( rrtm_play(0,nzt+8)   - rrtm_play(0,nzt+1)   )&
1234                              * ( rrtm_play(0,k) - rrtm_play(0,nzt+1) )
1235
1236             ENDDO
1237
1238!--          Linear interpolate to zw grid
1239             DO k = nzb+2, nzt+8
1240                rrtm_tlev(0,k)   = rrtm_tlay(0,k-1) + (rrtm_tlay(0,k) -        &
1241                                   rrtm_tlay(0,k-1))                           &
1242                                   / ( rrtm_play(0,k) - rrtm_play(0,k-1) )     &
1243                                   * ( rrtm_plev(0,k) - rrtm_play(0,k-1) )
1244             ENDDO
1245
1246
1247!
1248!--          Calculate liquid water path and cloud fraction for each column.
1249!--          Note that LWP is required in g/m² instead of kg/kg m.
1250             rrtm_cldfr  = 0.0_wp
1251             rrtm_reliq  = 0.0_wp
1252             rrtm_cliqwp = 0.0_wp
[1691]1253             rrtm_icld   = 0
[1585]1254
1255             DO k = nzb+1, nzt+1
[1691]1256                rrtm_cliqwp(0,k) =  ql(k,j,i) * 1000.0_wp *                    &
1257                                    (rrtm_plev(0,k) - rrtm_plev(0,k+1))        &
1258                                    * 100.0_wp / g 
[1585]1259
[1691]1260                IF ( rrtm_cliqwp(0,k) > 0.0_wp )  THEN
[1585]1261                   rrtm_cldfr(0,k) = 1.0_wp
[1691]1262                   IF ( rrtm_icld == 0 )  rrtm_icld = 1
[1585]1263
1264!
1265!--                Calculate cloud droplet effective radius
1266                   IF ( cloud_physics )  THEN
[1691]1267                      rrtm_reliq(0,k) = 1.0E6_wp * ( 3.0_wp * ql(k,j,i)        &
1268                                        * rho_surface                          &
1269                                        / ( 4.0_wp * pi * nc_const * rho_l )   &
1270                                        )**0.33333333333333_wp                 &
1271                                        * EXP( LOG( sigma_gc )**2 )
[1585]1272
1273                   ELSEIF ( cloud_droplets )  THEN
1274                      number_of_particles = prt_count(k,j,i)
1275
1276                      IF (number_of_particles <= 0)  CYCLE
1277                      particles => grid_particles(k,j,i)%particles(1:number_of_particles)
1278                      s_r2 = 0.0_wp
1279                      s_r3 = 0.0_wp
1280
1281                      DO  n = 1, number_of_particles
1282                         IF ( particles(n)%particle_mask )  THEN
1283                            s_r2 = s_r2 + particles(n)%radius**2 * &
1284                                   particles(n)%weight_factor
1285                            s_r3 = s_r3 + particles(n)%radius**3 * &
1286                                   particles(n)%weight_factor
1287                         ENDIF
1288                      ENDDO
1289
1290                      IF ( s_r2 > 0.0_wp )  rrtm_reliq(0,k) = s_r3 / s_r2
1291
1292                   ENDIF
1293
1294!
1295!--                Limit effective radius
[1691]1296                   IF ( rrtm_reliq(0,k) > 0.0_wp )  THEN
[1585]1297                      rrtm_reliq(0,k) = MAX(rrtm_reliq(0,k),2.5_wp)
1298                      rrtm_reliq(0,k) = MIN(rrtm_reliq(0,k),60.0_wp)
1299                  ENDIF
1300                ENDIF
1301             ENDDO
1302
1303!
1304!--          Set surface temperature
1305             rrtm_tsfc = pt(nzb,j,i) * (surface_pressure / 1000.0_wp )**0.286_wp
1306
1307             IF ( lw_radiation )  THEN
1308               CALL rrtmg_lw( 1, nzt_rad      , rrtm_icld    , rrtm_idrv      ,&
1309               rrtm_play       , rrtm_plev    , rrtm_tlay    , rrtm_tlev      ,&
1310               rrtm_tsfc       , rrtm_h2ovmr  , rrtm_o3vmr   , rrtm_co2vmr    ,&
1311               rrtm_ch4vmr     , rrtm_n2ovmr  , rrtm_o2vmr   , rrtm_cfc11vmr  ,&
1312               rrtm_cfc12vmr   , rrtm_cfc22vmr, rrtm_ccl4vmr , rrtm_emis      ,&
1313               rrtm_inflglw    , rrtm_iceflglw, rrtm_liqflglw, rrtm_cldfr     ,&
1314               rrtm_lw_taucld  , rrtm_cicewp  , rrtm_cliqwp  , rrtm_reice     ,& 
1315               rrtm_reliq      , rrtm_lw_tauaer,                               &
1316               rrtm_lwuflx     , rrtm_lwdflx  , rrtm_lwhr  ,                   &
[1691]1317               rrtm_lwuflxc    , rrtm_lwdflxc , rrtm_lwhrc ,                   &
1318               rrtm_lwuflx_dt  ,  rrtm_lwuflxc_dt )
[1585]1319
[1691]1320!
1321!--             Save fluxes
[1585]1322                DO k = nzb, nzt+1
1323                   rad_lw_in(k,j,i)  = rrtm_lwdflx(0,k)
1324                   rad_lw_out(k,j,i) = rrtm_lwuflx(0,k)
1325                ENDDO
1326
[1691]1327!
1328!--             Save heating rates (convert from K/d to K/h)
1329                DO k = nzb+1, nzt+1
1330                   rad_lw_hr(k,j,i)     = rrtm_lwhr(0,k)  * d_hours_day
1331                   rad_lw_cs_hr(k,j,i)  = rrtm_lwhrc(0,k) * d_hours_day
1332                ENDDO
[1585]1333
[1709]1334!
1335!--             Save change in LW heating rate
1336                rad_lw_out_change_0(j,i) = rrtm_lwuflx_dt(0,nzb)
1337
[1585]1338             ENDIF
1339
1340             IF ( sw_radiation .AND. sun_up )  THEN
1341                CALL rrtmg_sw( 1, nzt_rad      , rrtm_icld  , rrtm_iaer       ,&
1342               rrtm_play       , rrtm_plev    , rrtm_tlay  , rrtm_tlev        ,&
1343               rrtm_tsfc       , rrtm_h2ovmr  , rrtm_o3vmr , rrtm_co2vmr      ,&
1344               rrtm_ch4vmr     , rrtm_n2ovmr  , rrtm_o2vmr , rrtm_asdir(:,j,i),&
1345               rrtm_asdif(:,j,i), rrtm_aldir(:,j,i), rrtm_aldif(:,j,i), zenith,&
1346               0.0_wp          , day          , solar_constant,   rrtm_inflgsw,&
1347               rrtm_iceflgsw   , rrtm_liqflgsw, rrtm_cldfr , rrtm_sw_taucld   ,&
1348               rrtm_sw_ssacld  , rrtm_sw_asmcld, rrtm_sw_fsfcld, rrtm_cicewp  ,&
1349               rrtm_cliqwp     , rrtm_reice   , rrtm_reliq , rrtm_sw_tauaer   ,&
1350               rrtm_sw_ssaaer     , rrtm_sw_asmaer  , rrtm_sw_ecaer ,          &
1351               rrtm_swuflx     , rrtm_swdflx  , rrtm_swhr  ,                   &
1352               rrtm_swuflxc    , rrtm_swdflxc , rrtm_swhrc )
1353 
[1691]1354!
1355!--             Save fluxes
[1585]1356                DO k = nzb, nzt+1
1357                   rad_sw_in(k,j,i)  = rrtm_swdflx(0,k)
1358                   rad_sw_out(k,j,i) = rrtm_swuflx(0,k)
1359                ENDDO
[1691]1360
1361!
1362!--             Save heating rates (convert from K/d to K/s)
1363                DO k = nzb+1, nzt+1
1364                   rad_sw_hr(k,j,i)     = rrtm_swhr(0,k)  * d_hours_day
1365                   rad_sw_cs_hr(k,j,i)  = rrtm_swhrc(0,k) * d_hours_day
1366                ENDDO
1367
[1585]1368             ENDIF
1369
1370!
1371!--          Calculate surface net radiation
1372             rad_net(j,i) = rad_sw_in(nzb,j,i) - rad_sw_out(nzb,j,i)           &
1373                            + rad_lw_in(nzb,j,i) - rad_lw_out(nzb,j,i)
1374
1375          ENDDO
1376       ENDDO
1377
1378       CALL exchange_horiz( rad_lw_in,  nbgp )
1379       CALL exchange_horiz( rad_lw_out, nbgp )
[1691]1380       CALL exchange_horiz( rad_lw_hr,    nbgp )
1381       CALL exchange_horiz( rad_lw_cs_hr, nbgp )
1382
[1585]1383       CALL exchange_horiz( rad_sw_in,  nbgp )
1384       CALL exchange_horiz( rad_sw_out, nbgp ) 
[1691]1385       CALL exchange_horiz( rad_sw_hr,    nbgp )
1386       CALL exchange_horiz( rad_sw_cs_hr, nbgp )
1387
[1585]1388       CALL exchange_horiz_2d( rad_net, nbgp )
[1709]1389       CALL exchange_horiz_2d( rad_lw_out_change_0, nbgp )
[1585]1390#endif
1391
1392    END SUBROUTINE radiation_rrtmg
1393
1394
1395!------------------------------------------------------------------------------!
1396! Description:
1397! ------------
[1682]1398!> Calculate the cosine of the zenith angle (variable is called zenith)
[1585]1399!------------------------------------------------------------------------------!
1400    SUBROUTINE calc_zenith
1401
1402       IMPLICIT NONE
1403
[1682]1404       REAL(wp) ::  declination,  & !< solar declination angle
1405                    hour_angle      !< solar hour angle
[1585]1406!
[1496]1407!--    Calculate current day and time based on the initial values and simulation
1408!--    time
[1585]1409       day = day_init + INT(FLOOR( (time_utc_init + time_since_reference_point)    &
1410                               / 86400.0_wp ), KIND=iwp)
[1496]1411       time_utc = MOD((time_utc_init + time_since_reference_point), 86400.0_wp)
1412
1413
1414!
1415!--    Calculate solar declination and hour angle   
[1585]1416       declination = ASIN( decl_1 * SIN(decl_2 * REAL(day, KIND=wp) - decl_3) )
[1496]1417       hour_angle  = 2.0_wp * pi * (time_utc / 86400.0_wp) + lon - pi
1418
1419!
1420!--    Calculate zenith angle
[1585]1421       zenith(0) = SIN(lat) * SIN(declination) + COS(lat) * COS(declination)      &
[1496]1422                                            * COS(hour_angle)
[1585]1423       zenith(0) = MAX(0.0_wp,zenith(0))
[1496]1424
1425!
[1585]1426!--    Check if the sun is up (otheriwse shortwave calculations can be skipped)
[1691]1427       IF ( zenith(0) > 0.0_wp )  THEN
[1585]1428          sun_up = .TRUE.
1429       ELSE
1430          sun_up = .FALSE.
1431       END IF
[1496]1432
[1585]1433    END SUBROUTINE calc_zenith
1434
[1606]1435#if defined ( __rrtmg ) && defined ( __netcdf )
[1585]1436!------------------------------------------------------------------------------!
1437! Description:
1438! ------------
[1682]1439!> Calculates surface albedo components based on Briegleb (1992) and
1440!> Briegleb et al. (1986)
[1585]1441!------------------------------------------------------------------------------!
1442    SUBROUTINE calc_albedo
1443
1444        IMPLICIT NONE
1445
1446        IF ( sun_up )  THEN
[1496]1447!
[1585]1448!--        Ocean
1449           IF ( albedo_type == 1 )  THEN
1450              rrtm_aldir(0,:,:) = 0.026_wp / ( zenith(0)**1.7_wp + 0.065_wp )  &
1451                                  + 0.15_wp * ( zenith(0) - 0.1_wp )           &
1452                                            * ( zenith(0) - 0.5_wp )           &
1453                                            * ( zenith(0) - 1.0_wp )
1454              rrtm_asdir(0,:,:) = rrtm_aldir(0,:,:)
1455!
1456!--        Snow
1457           ELSEIF ( albedo_type == 16 )  THEN
1458              IF ( zenith(0) < 0.5_wp )  THEN
1459                 rrtm_aldir(0,:,:) = 0.5_wp * (1.0_wp - aldif)                 &
1460                                     * ( 3.0_wp / (1.0_wp + 4.0_wp             &
1461                                     * zenith(0))) - 1.0_wp
1462                 rrtm_asdir(0,:,:) = 0.5_wp * (1.0_wp - asdif)                 &
1463                                     * ( 3.0_wp / (1.0_wp + 4.0_wp             &
1464                                     * zenith(0))) - 1.0_wp
[1496]1465
[1585]1466                 rrtm_aldir(0,:,:) = MIN(0.98_wp, rrtm_aldir(0,:,:))
1467                 rrtm_asdir(0,:,:) = MIN(0.98_wp, rrtm_asdir(0,:,:))
1468              ELSE
1469                 rrtm_aldir(0,:,:) = aldif
1470                 rrtm_asdir(0,:,:) = asdif
1471              ENDIF
[1496]1472!
[1585]1473!--        Sea ice
1474           ELSEIF ( albedo_type == 15 )  THEN
1475                 rrtm_aldir(0,:,:) = aldif
1476                 rrtm_asdir(0,:,:) = asdif
[1788]1477
[1585]1478!
[1788]1479!--        Asphalt
1480           ELSEIF ( albedo_type == 17 )  THEN
1481                 rrtm_aldir(0,:,:) = aldif
1482                 rrtm_asdir(0,:,:) = asdif
1483!
[1585]1484!--        Land surfaces
1485           ELSE
1486              SELECT CASE ( albedo_type )
[1496]1487
[1585]1488!
1489!--              Surface types with strong zenith dependence
1490                 CASE ( 1, 2, 3, 4, 11, 12, 13 )
1491                    rrtm_aldir(0,:,:) = aldif * 1.4_wp /                       &
1492                                        (1.0_wp + 0.8_wp * zenith(0))
1493                    rrtm_asdir(0,:,:) = asdif * 1.4_wp /                       &
1494                                        (1.0_wp + 0.8_wp * zenith(0))
1495!
1496!--              Surface types with weak zenith dependence
1497                 CASE ( 5, 6, 7, 8, 9, 10, 14 )
1498                    rrtm_aldir(0,:,:) = aldif * 1.1_wp /                       &
1499                                        (1.0_wp + 0.2_wp * zenith(0))
1500                    rrtm_asdir(0,:,:) = asdif * 1.1_wp /                       &
1501                                        (1.0_wp + 0.2_wp * zenith(0))
[1496]1502
[1585]1503                 CASE DEFAULT
1504
1505              END SELECT
1506           ENDIF
1507!
1508!--        Diffusive albedo is taken from Table 2
1509           rrtm_aldif(0,:,:) = aldif
1510           rrtm_asdif(0,:,:) = asdif
1511
1512        ELSE
1513
1514           rrtm_aldir(0,:,:) = 0.0_wp
1515           rrtm_asdir(0,:,:) = 0.0_wp
1516           rrtm_aldif(0,:,:) = 0.0_wp
1517           rrtm_asdif(0,:,:) = 0.0_wp
1518        ENDIF
1519    END SUBROUTINE calc_albedo
1520
1521!------------------------------------------------------------------------------!
1522! Description:
1523! ------------
[1682]1524!> Read sounding data (pressure and temperature) from RADIATION_DATA.
[1585]1525!------------------------------------------------------------------------------!
1526    SUBROUTINE read_sounding_data
1527
1528       IMPLICIT NONE
1529
[1691]1530       INTEGER(iwp) :: id,           & !< NetCDF id of input file
1531                       id_dim_zrad,  & !< pressure level id in the NetCDF file
1532                       id_var,       & !< NetCDF variable id
1533                       k,            & !< loop index
1534                       nz_snd,       & !< number of vertical levels in the sounding data
1535                       nz_snd_start, & !< start vertical index for sounding data to be used
1536                       nz_snd_end      !< end vertical index for souding data to be used
[1585]1537
[1691]1538       REAL(wp) :: t_surface           !< actual surface temperature
[1585]1539
[1691]1540       REAL(wp), DIMENSION(:), ALLOCATABLE ::  hyp_snd_tmp, & !< temporary hydrostatic pressure profile (sounding)
1541                                               t_snd_tmp      !< temporary temperature profile (sounding)
[1585]1542
1543!
1544!--    In case of updates, deallocate arrays first (sufficient to check one
1545!--    array as the others are automatically allocated). This is required
1546!--    because nzt_rad might change during the update
1547       IF ( ALLOCATED ( hyp_snd ) )  THEN
1548          DEALLOCATE( hyp_snd )
1549          DEALLOCATE( t_snd )
1550          DEALLOCATE( q_snd  )
1551          DEALLOCATE ( rrtm_play )
1552          DEALLOCATE ( rrtm_plev )
1553          DEALLOCATE ( rrtm_tlay )
1554          DEALLOCATE ( rrtm_tlev )
[1691]1555
[1585]1556          DEALLOCATE ( rrtm_h2ovmr )
1557          DEALLOCATE ( rrtm_cicewp )
1558          DEALLOCATE ( rrtm_cldfr )
1559          DEALLOCATE ( rrtm_cliqwp )
1560          DEALLOCATE ( rrtm_reice )
1561          DEALLOCATE ( rrtm_reliq )
1562          DEALLOCATE ( rrtm_lw_taucld )
1563          DEALLOCATE ( rrtm_lw_tauaer )
[1691]1564
[1585]1565          DEALLOCATE ( rrtm_lwdflx  )
[1691]1566          DEALLOCATE ( rrtm_lwdflxc )
[1585]1567          DEALLOCATE ( rrtm_lwuflx  )
[1691]1568          DEALLOCATE ( rrtm_lwuflxc )
1569          DEALLOCATE ( rrtm_lwuflx_dt )
1570          DEALLOCATE ( rrtm_lwuflxc_dt )
[1585]1571          DEALLOCATE ( rrtm_lwhr  )
1572          DEALLOCATE ( rrtm_lwhrc )
[1691]1573
[1585]1574          DEALLOCATE ( rrtm_sw_taucld )
1575          DEALLOCATE ( rrtm_sw_ssacld )
1576          DEALLOCATE ( rrtm_sw_asmcld )
1577          DEALLOCATE ( rrtm_sw_fsfcld )
1578          DEALLOCATE ( rrtm_sw_tauaer )
1579          DEALLOCATE ( rrtm_sw_ssaaer )
1580          DEALLOCATE ( rrtm_sw_asmaer ) 
[1691]1581          DEALLOCATE ( rrtm_sw_ecaer )   
1582 
[1585]1583          DEALLOCATE ( rrtm_swdflx  )
[1691]1584          DEALLOCATE ( rrtm_swdflxc )
[1585]1585          DEALLOCATE ( rrtm_swuflx  )
[1691]1586          DEALLOCATE ( rrtm_swuflxc )
[1585]1587          DEALLOCATE ( rrtm_swhr  )
1588          DEALLOCATE ( rrtm_swhrc )
[1691]1589
[1585]1590       ENDIF
1591
1592!
1593!--    Open file for reading
1594       nc_stat = NF90_OPEN( rrtm_input_file, NF90_NOWRITE, id )
[1783]1595       CALL netcdf_handle_error_rad( 'read_sounding_data', 549 )
[1585]1596
1597!
1598!--    Inquire dimension of z axis and save in nz_snd
1599       nc_stat = NF90_INQ_DIMID( id, "Pressure", id_dim_zrad )
1600       nc_stat = NF90_INQUIRE_DIMENSION( id, id_dim_zrad, len = nz_snd )
[1783]1601       CALL netcdf_handle_error_rad( 'read_sounding_data', 551 )
[1585]1602
1603!
1604! !--    Allocate temporary array for storing pressure data
[1701]1605       ALLOCATE( hyp_snd_tmp(1:nz_snd) )
[1585]1606       hyp_snd_tmp = 0.0_wp
1607
1608
1609!--    Read pressure from file
1610       nc_stat = NF90_INQ_VARID( id, "Pressure", id_var )
[1691]1611       nc_stat = NF90_GET_VAR( id, id_var, hyp_snd_tmp(:), start = (/1/),      &
[1585]1612                               count = (/nz_snd/) )
[1783]1613       CALL netcdf_handle_error_rad( 'read_sounding_data', 552 )
[1585]1614
1615!
1616!--    Allocate temporary array for storing temperature data
[1701]1617       ALLOCATE( t_snd_tmp(1:nz_snd) )
[1585]1618       t_snd_tmp = 0.0_wp
1619
1620!
1621!--    Read temperature from file
1622       nc_stat = NF90_INQ_VARID( id, "ReferenceTemperature", id_var )
[1691]1623       nc_stat = NF90_GET_VAR( id, id_var, t_snd_tmp(:), start = (/1/),        &
[1585]1624                               count = (/nz_snd/) )
[1783]1625       CALL netcdf_handle_error_rad( 'read_sounding_data', 553 )
[1585]1626
1627!
1628!--    Calculate start of sounding data
1629       nz_snd_start = nz_snd + 1
[1701]1630       nz_snd_end   = nz_snd + 1
[1585]1631
1632!
1633!--    Start filling vertical dimension at 10hPa above the model domain (hyp is
1634!--    in Pa, hyp_snd in hPa).
1635       DO  k = 1, nz_snd
[1691]1636          IF ( hyp_snd_tmp(k) < ( hyp(nzt+1) - 1000.0_wp) * 0.01_wp )  THEN
[1585]1637             nz_snd_start = k
1638             EXIT
1639          END IF
1640       END DO
1641
[1691]1642       IF ( nz_snd_start <= nz_snd )  THEN
[1701]1643          nz_snd_end = nz_snd
[1585]1644       END IF
1645
1646
1647!
1648!--    Calculate of total grid points for RRTMG calculations
[1701]1649       nzt_rad = nzt + nz_snd_end - nz_snd_start + 1
[1585]1650
1651!
1652!--    Save data above LES domain in hyp_snd, t_snd and q_snd
1653!--    Note: q_snd_tmp is not calculated at the moment (dry residual atmosphere)
1654       ALLOCATE( hyp_snd(nzb+1:nzt_rad) )
1655       ALLOCATE( t_snd(nzb+1:nzt_rad)   )
1656       ALLOCATE( q_snd(nzb+1:nzt_rad)   )
1657       hyp_snd = 0.0_wp
1658       t_snd = 0.0_wp
1659       q_snd = 0.0_wp
1660
[1757]1661       hyp_snd(nzt+2:nzt_rad) = hyp_snd_tmp(nz_snd_start+1:nz_snd_end)
1662       t_snd(nzt+2:nzt_rad)   = t_snd_tmp(nz_snd_start+1:nz_snd_end)
[1585]1663
1664       nc_stat = NF90_CLOSE( id )
1665
1666!
1667!--    Calculate pressure levels on zu and zw grid. Sounding data is added at
1668!--    top of the LES domain. This routine does not consider horizontal or
1669!--    vertical variability of pressure and temperature
1670       ALLOCATE ( rrtm_play(0:0,nzb+1:nzt_rad+1)   )
1671       ALLOCATE ( rrtm_plev(0:0,nzb+1:nzt_rad+2)   )
1672
[1691]1673       t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**0.286_wp
[1585]1674       DO k = nzb+1, nzt+1
1675          rrtm_play(0,k) = hyp(k) * 0.01_wp
1676          rrtm_plev(0,k) = surface_pressure * ( (t_surface - g/cp * zw(k-1)) / &
1677                         t_surface )**(1.0_wp/0.286_wp)
1678       ENDDO
1679
1680       DO k = nzt+2, nzt_rad
1681          rrtm_play(0,k) = hyp_snd(k)
1682          rrtm_plev(0,k) = 0.5_wp * ( rrtm_play(0,k) + rrtm_play(0,k-1) )
1683       ENDDO
1684       rrtm_plev(0,nzt_rad+1) = MAX( 0.5 * hyp_snd(nzt_rad),                   &
1685                                   1.5 * hyp_snd(nzt_rad)                      &
1686                                 - 0.5 * hyp_snd(nzt_rad-1) )
1687       rrtm_plev(0,nzt_rad+2)  = MIN( 1.0E-4_wp,                               &
1688                                      0.25_wp * rrtm_plev(0,nzt_rad+1) )
1689
1690       rrtm_play(0,nzt_rad+1) = 0.5 * rrtm_plev(0,nzt_rad+1)
1691
1692!
1693!--    Calculate temperature/humidity levels at top of the LES domain.
1694!--    Currently, the temperature is taken from sounding data (might lead to a
1695!--    temperature jump at interface. To do: Humidity is currently not
1696!--    calculated above the LES domain.
1697       ALLOCATE ( rrtm_tlay(0:0,nzb+1:nzt_rad+1)   )
1698       ALLOCATE ( rrtm_tlev(0:0,nzb+1:nzt_rad+2)   )
1699       ALLOCATE ( rrtm_h2ovmr(0:0,nzb+1:nzt_rad+1) )
1700
1701       DO k = nzt+8, nzt_rad
1702          rrtm_tlay(0,k)   = t_snd(k)
1703          rrtm_h2ovmr(0,k) = q_snd(k)
1704       ENDDO
[1691]1705       rrtm_tlay(0,nzt_rad+1) = 2.0_wp * rrtm_tlay(0,nzt_rad)                 &
1706                                - rrtm_tlay(0,nzt_rad-1)
[1585]1707       DO k = nzt+9, nzt_rad+1
1708          rrtm_tlev(0,k)   = rrtm_tlay(0,k-1) + (rrtm_tlay(0,k)                &
1709                             - rrtm_tlay(0,k-1))                               &
1710                             / ( rrtm_play(0,k) - rrtm_play(0,k-1) )           &
1711                             * ( rrtm_plev(0,k) - rrtm_play(0,k-1) )
1712       ENDDO
1713       rrtm_h2ovmr(0,nzt_rad+1) = rrtm_h2ovmr(0,nzt_rad)
1714
1715       rrtm_tlev(0,nzt_rad+2)   = 2.0_wp * rrtm_tlay(0,nzt_rad+1)              &
1716                                  - rrtm_tlev(0,nzt_rad)
1717!
1718!--    Allocate remaining RRTMG arrays
1719       ALLOCATE ( rrtm_cicewp(0:0,nzb+1:nzt_rad+1) )
1720       ALLOCATE ( rrtm_cldfr(0:0,nzb+1:nzt_rad+1) )
1721       ALLOCATE ( rrtm_cliqwp(0:0,nzb+1:nzt_rad+1) )
1722       ALLOCATE ( rrtm_reice(0:0,nzb+1:nzt_rad+1) )
1723       ALLOCATE ( rrtm_reliq(0:0,nzb+1:nzt_rad+1) )
1724       ALLOCATE ( rrtm_lw_taucld(1:nbndlw+1,0:0,nzb+1:nzt_rad+1) )
1725       ALLOCATE ( rrtm_lw_tauaer(0:0,nzb+1:nzt_rad+1,1:nbndlw+1) )
1726       ALLOCATE ( rrtm_sw_taucld(1:nbndsw+1,0:0,nzb+1:nzt_rad+1) )
1727       ALLOCATE ( rrtm_sw_ssacld(1:nbndsw+1,0:0,nzb+1:nzt_rad+1) )
1728       ALLOCATE ( rrtm_sw_asmcld(1:nbndsw+1,0:0,nzb+1:nzt_rad+1) )
1729       ALLOCATE ( rrtm_sw_fsfcld(1:nbndsw+1,0:0,nzb+1:nzt_rad+1) )
1730       ALLOCATE ( rrtm_sw_tauaer(0:0,nzb+1:nzt_rad+1,1:nbndsw+1) )
1731       ALLOCATE ( rrtm_sw_ssaaer(0:0,nzb+1:nzt_rad+1,1:nbndsw+1) )
1732       ALLOCATE ( rrtm_sw_asmaer(0:0,nzb+1:nzt_rad+1,1:nbndsw+1) ) 
1733       ALLOCATE ( rrtm_sw_ecaer(0:0,nzb+1:nzt_rad+1,1:naerec+1) )   
1734
1735!
1736!--    The ice phase is currently not considered in PALM
1737       rrtm_cicewp = 0.0_wp
1738       rrtm_reice  = 0.0_wp
1739
1740!
1741!--    Set other parameters (move to NAMELIST parameters in the future)
1742       rrtm_lw_tauaer = 0.0_wp
1743       rrtm_lw_taucld = 0.0_wp
1744       rrtm_sw_taucld = 0.0_wp
1745       rrtm_sw_ssacld = 0.0_wp
1746       rrtm_sw_asmcld = 0.0_wp
1747       rrtm_sw_fsfcld = 0.0_wp
1748       rrtm_sw_tauaer = 0.0_wp
1749       rrtm_sw_ssaaer = 0.0_wp
1750       rrtm_sw_asmaer = 0.0_wp
1751       rrtm_sw_ecaer  = 0.0_wp
1752
1753
1754       ALLOCATE ( rrtm_swdflx(0:0,nzb:nzt_rad+1)  )
1755       ALLOCATE ( rrtm_swuflx(0:0,nzb:nzt_rad+1)  )
1756       ALLOCATE ( rrtm_swhr(0:0,nzb+1:nzt_rad+1)  )
1757       ALLOCATE ( rrtm_swuflxc(0:0,nzb:nzt_rad+1) )
1758       ALLOCATE ( rrtm_swdflxc(0:0,nzb:nzt_rad+1) )
1759       ALLOCATE ( rrtm_swhrc(0:0,nzb+1:nzt_rad+1) )
1760
1761       rrtm_swdflx  = 0.0_wp
1762       rrtm_swuflx  = 0.0_wp
1763       rrtm_swhr    = 0.0_wp 
1764       rrtm_swuflxc = 0.0_wp
1765       rrtm_swdflxc = 0.0_wp
1766       rrtm_swhrc   = 0.0_wp
1767
1768       ALLOCATE ( rrtm_lwdflx(0:0,nzb:nzt_rad+1)  )
1769       ALLOCATE ( rrtm_lwuflx(0:0,nzb:nzt_rad+1)  )
1770       ALLOCATE ( rrtm_lwhr(0:0,nzb+1:nzt_rad+1)  )
1771       ALLOCATE ( rrtm_lwuflxc(0:0,nzb:nzt_rad+1) )
1772       ALLOCATE ( rrtm_lwdflxc(0:0,nzb:nzt_rad+1) )
1773       ALLOCATE ( rrtm_lwhrc(0:0,nzb+1:nzt_rad+1) )
1774
1775       rrtm_lwdflx  = 0.0_wp
1776       rrtm_lwuflx  = 0.0_wp
1777       rrtm_lwhr    = 0.0_wp 
1778       rrtm_lwuflxc = 0.0_wp
1779       rrtm_lwdflxc = 0.0_wp
1780       rrtm_lwhrc   = 0.0_wp
1781
[1691]1782       ALLOCATE ( rrtm_lwuflx_dt(0:0,nzb:nzt_rad+1) )
1783       ALLOCATE ( rrtm_lwuflxc_dt(0:0,nzb:nzt_rad+1) )
[1585]1784
[1709]1785       rrtm_lwuflx_dt = 0.0_wp
[1691]1786       rrtm_lwuflxc_dt = 0.0_wp
1787
[1585]1788    END SUBROUTINE read_sounding_data
1789
1790
1791!------------------------------------------------------------------------------!
1792! Description:
1793! ------------
[1682]1794!> Read trace gas data from file
[1585]1795!------------------------------------------------------------------------------!
1796    SUBROUTINE read_trace_gas_data
1797
1798       USE rrsw_ncpar
1799
1800       IMPLICIT NONE
1801
[1691]1802       INTEGER(iwp), PARAMETER :: num_trace_gases = 9 !< number of trace gases (absorbers)
[1585]1803
[1691]1804       CHARACTER(LEN=5), DIMENSION(num_trace_gases), PARAMETER ::              & !< trace gas names
[1585]1805           trace_names = (/'O3   ', 'CO2  ', 'CH4  ', 'N2O  ', 'O2   ',        &
1806                           'CFC11', 'CFC12', 'CFC22', 'CCL4 '/)
1807
[1691]1808       INTEGER(iwp) :: id,     & !< NetCDF id
1809                       k,      & !< loop index
1810                       m,      & !< loop index
1811                       n,      & !< loop index
1812                       nabs,   & !< number of absorbers
1813                       np,     & !< number of pressure levels
1814                       id_abs, & !< NetCDF id of the respective absorber
1815                       id_dim, & !< NetCDF id of asborber's dimension
1816                       id_var    !< NetCDf id ot the absorber
[1585]1817
1818       REAL(wp) :: p_mls_l, p_mls_u, p_wgt_l, p_wgt_u, p_mls_m
1819
1820
[1682]1821       REAL(wp), DIMENSION(:), ALLOCATABLE   ::  p_mls,         & !< pressure levels for the absorbers
1822                                                 rrtm_play_tmp, & !< temporary array for pressure zu-levels
1823                                                 rrtm_plev_tmp, & !< temporary array for pressure zw-levels
1824                                                 trace_path_tmp   !< temporary array for storing trace gas path data
[1585]1825
[1682]1826       REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  trace_mls,      & !< array for storing the absorber amounts
1827                                                 trace_mls_path, & !< array for storing trace gas path data
1828                                                 trace_mls_tmp     !< temporary array for storing trace gas data
[1585]1829
1830
1831!
1832!--    In case of updates, deallocate arrays first (sufficient to check one
1833!--    array as the others are automatically allocated)
1834       IF ( ALLOCATED ( rrtm_o3vmr ) )  THEN
1835          DEALLOCATE ( rrtm_o3vmr  )
1836          DEALLOCATE ( rrtm_co2vmr )
1837          DEALLOCATE ( rrtm_ch4vmr )
1838          DEALLOCATE ( rrtm_n2ovmr )
1839          DEALLOCATE ( rrtm_o2vmr  )
1840          DEALLOCATE ( rrtm_cfc11vmr )
1841          DEALLOCATE ( rrtm_cfc12vmr )
1842          DEALLOCATE ( rrtm_cfc22vmr )
1843          DEALLOCATE ( rrtm_ccl4vmr  )
1844       ENDIF
1845
1846!
1847!--    Allocate trace gas profiles
1848       ALLOCATE ( rrtm_o3vmr(0:0,1:nzt_rad+1)  )
1849       ALLOCATE ( rrtm_co2vmr(0:0,1:nzt_rad+1) )
1850       ALLOCATE ( rrtm_ch4vmr(0:0,1:nzt_rad+1) )
1851       ALLOCATE ( rrtm_n2ovmr(0:0,1:nzt_rad+1) )
1852       ALLOCATE ( rrtm_o2vmr(0:0,1:nzt_rad+1)  )
1853       ALLOCATE ( rrtm_cfc11vmr(0:0,1:nzt_rad+1) )
1854       ALLOCATE ( rrtm_cfc12vmr(0:0,1:nzt_rad+1) )
1855       ALLOCATE ( rrtm_cfc22vmr(0:0,1:nzt_rad+1) )
1856       ALLOCATE ( rrtm_ccl4vmr(0:0,1:nzt_rad+1)  )
1857
1858!
1859!--    Open file for reading
1860       nc_stat = NF90_OPEN( rrtm_input_file, NF90_NOWRITE, id )
[1783]1861       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 549 )
[1585]1862!
1863!--    Inquire dimension ids and dimensions
1864       nc_stat = NF90_INQ_DIMID( id, "Pressure", id_dim )
[1783]1865       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1866       nc_stat = NF90_INQUIRE_DIMENSION( id, id_dim, len = np) 
[1783]1867       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1868
1869       nc_stat = NF90_INQ_DIMID( id, "Absorber", id_dim )
[1783]1870       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1871       nc_stat = NF90_INQUIRE_DIMENSION( id, id_dim, len = nabs ) 
[1783]1872       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1873   
1874
1875!
1876!--    Allocate pressure, and trace gas arrays     
1877       ALLOCATE( p_mls(1:np) )
1878       ALLOCATE( trace_mls(1:num_trace_gases,1:np) ) 
1879       ALLOCATE( trace_mls_tmp(1:nabs,1:np) ) 
1880
1881
1882       nc_stat = NF90_INQ_VARID( id, "Pressure", id_var )
[1783]1883       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1884       nc_stat = NF90_GET_VAR( id, id_var, p_mls )
[1783]1885       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1886
1887       nc_stat = NF90_INQ_VARID( id, "AbsorberAmountMLS", id_var )
[1783]1888       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1889       nc_stat = NF90_GET_VAR( id, id_var, trace_mls_tmp )
[1783]1890       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 550 )
[1585]1891
1892
1893!
1894!--    Write absorber amounts (mls) to trace_mls
1895       DO n = 1, num_trace_gases
1896          CALL getAbsorberIndex( TRIM( trace_names(n) ), id_abs )
1897
1898          trace_mls(n,1:np) = trace_mls_tmp(id_abs,1:np)
1899
1900!
1901!--       Replace missing values by zero
1902          WHERE ( trace_mls(n,:) > 2.0_wp ) 
1903             trace_mls(n,:) = 0.0_wp
1904          END WHERE
1905       END DO
1906
1907       DEALLOCATE ( trace_mls_tmp )
1908
1909       nc_stat = NF90_CLOSE( id )
[1783]1910       CALL netcdf_handle_error_rad( 'read_trace_gas_data', 551 )
[1585]1911
1912!
1913!--    Add extra pressure level for calculations of the trace gas paths
1914       ALLOCATE ( rrtm_play_tmp(1:nzt_rad+1) )
1915       ALLOCATE ( rrtm_plev_tmp(1:nzt_rad+2) )
1916
1917       rrtm_play_tmp(1:nzt_rad)   = rrtm_play(0,1:nzt_rad) 
1918       rrtm_plev_tmp(1:nzt_rad+1) = rrtm_plev(0,1:nzt_rad+1)
1919       rrtm_play_tmp(nzt_rad+1)   = rrtm_plev(0,nzt_rad+1) * 0.5_wp
1920       rrtm_plev_tmp(nzt_rad+2)   = MIN( 1.0E-4_wp, 0.25_wp                    &
1921                                         * rrtm_plev(0,nzt_rad+1) )
1922 
1923!
1924!--    Calculate trace gas path (zero at surface) with interpolation to the
1925!--    sounding levels
1926       ALLOCATE ( trace_mls_path(1:nzt_rad+2,1:num_trace_gases) )
1927
1928       trace_mls_path(nzb+1,:) = 0.0_wp
1929       
1930       DO k = nzb+2, nzt_rad+2
1931          DO m = 1, num_trace_gases
1932             trace_mls_path(k,m) = trace_mls_path(k-1,m)
1933
1934!
1935!--          When the pressure level is higher than the trace gas pressure
1936!--          level, assume that
[1691]1937             IF ( rrtm_plev_tmp(k-1) > p_mls(1) )  THEN             
[1585]1938               
1939                trace_mls_path(k,m) = trace_mls_path(k,m) + trace_mls(m,1)     &
1940                                      * ( rrtm_plev_tmp(k-1)                   &
1941                                          - MAX( p_mls(1), rrtm_plev_tmp(k) )  &
1942                                        ) / g
1943             ENDIF
1944
1945!
1946!--          Integrate for each sounding level from the contributing p_mls
1947!--          levels
1948             DO n = 2, np
1949!
1950!--             Limit p_mls so that it is within the model level
1951                p_mls_u = MIN( rrtm_plev_tmp(k-1),                             &
1952                          MAX( rrtm_plev_tmp(k), p_mls(n) ) )
1953                p_mls_l = MIN( rrtm_plev_tmp(k-1),                             &
1954                          MAX( rrtm_plev_tmp(k), p_mls(n-1) ) )
1955
[1691]1956                IF ( p_mls_l > p_mls_u )  THEN
[1585]1957
1958!
1959!--                Calculate weights for interpolation
1960                   p_mls_m = 0.5_wp * (p_mls_l + p_mls_u)
1961                   p_wgt_u = (p_mls(n-1) - p_mls_m) / (p_mls(n-1) - p_mls(n))
1962                   p_wgt_l = (p_mls_m - p_mls(n))   / (p_mls(n-1) - p_mls(n))
1963
1964!
1965!--                Add level to trace gas path
1966                   trace_mls_path(k,m) = trace_mls_path(k,m)                   &
1967                                         +  ( p_wgt_u * trace_mls(m,n)         &
1968                                            + p_wgt_l * trace_mls(m,n-1) )     &
[1691]1969                                         * (p_mls_l - p_mls_u) / g
[1585]1970                ENDIF
1971             ENDDO
1972
[1691]1973             IF ( rrtm_plev_tmp(k) < p_mls(np) )  THEN
[1585]1974                trace_mls_path(k,m) = trace_mls_path(k,m) + trace_mls(m,np)    &
1975                                      * ( MIN( rrtm_plev_tmp(k-1), p_mls(np) ) &
1976                                          - rrtm_plev_tmp(k)                   &
1977                                        ) / g 
1978             ENDIF 
[1496]1979          ENDDO
1980       ENDDO
1981
1982
[1585]1983!
1984!--    Prepare trace gas path profiles
1985       ALLOCATE ( trace_path_tmp(1:nzt_rad+1) )
[1496]1986
[1585]1987       DO m = 1, num_trace_gases
1988
1989          trace_path_tmp(1:nzt_rad+1) = ( trace_mls_path(2:nzt_rad+2,m)        &
1990                                       - trace_mls_path(1:nzt_rad+1,m) ) * g   &
1991                                       / ( rrtm_plev_tmp(1:nzt_rad+1)          &
1992                                       - rrtm_plev_tmp(2:nzt_rad+2) )
1993
1994!
1995!--       Save trace gas paths to the respective arrays
1996          SELECT CASE ( TRIM( trace_names(m) ) )
1997
1998             CASE ( 'O3' )
1999
2000                rrtm_o3vmr(0,:) = trace_path_tmp(:)
2001
2002             CASE ( 'CO2' )
2003
2004                rrtm_co2vmr(0,:) = trace_path_tmp(:)
2005
2006             CASE ( 'CH4' )
2007
2008                rrtm_ch4vmr(0,:) = trace_path_tmp(:)
2009
2010             CASE ( 'N2O' )
2011
2012                rrtm_n2ovmr(0,:) = trace_path_tmp(:)
2013
2014             CASE ( 'O2' )
2015
2016                rrtm_o2vmr(0,:) = trace_path_tmp(:)
2017
2018             CASE ( 'CFC11' )
2019
2020                rrtm_cfc11vmr(0,:) = trace_path_tmp(:)
2021
2022             CASE ( 'CFC12' )
2023
2024                rrtm_cfc12vmr(0,:) = trace_path_tmp(:)
2025
2026             CASE ( 'CFC22' )
2027
2028                rrtm_cfc22vmr(0,:) = trace_path_tmp(:)
2029
2030             CASE ( 'CCL4' )
2031
2032                rrtm_ccl4vmr(0,:) = trace_path_tmp(:)
2033
2034             CASE DEFAULT
2035
2036          END SELECT
2037
2038       ENDDO
2039
2040       DEALLOCATE ( trace_path_tmp )
2041       DEALLOCATE ( trace_mls_path )
2042       DEALLOCATE ( rrtm_play_tmp )
2043       DEALLOCATE ( rrtm_plev_tmp )
2044       DEALLOCATE ( trace_mls )
2045       DEALLOCATE ( p_mls )
2046
2047    END SUBROUTINE read_trace_gas_data
2048
[1826]2049
[1783]2050    SUBROUTINE netcdf_handle_error_rad( routine_name, errno )
2051
2052       USE control_parameters,                                                 &
2053           ONLY:  message_string
2054
2055       USE NETCDF
2056
2057       USE pegrid
2058
2059       IMPLICIT NONE
2060
2061       CHARACTER(LEN=6) ::  message_identifier
2062       CHARACTER(LEN=*) ::  routine_name
2063
2064       INTEGER(iwp) ::  errno
2065
2066       IF ( nc_stat /= NF90_NOERR )  THEN
2067
2068          WRITE( message_identifier, '(''NC'',I4.4)' )  errno
2069          message_string = TRIM( NF90_STRERROR( nc_stat ) )
2070
2071          CALL message( routine_name, message_identifier, 2, 2, 0, 6, 1 )
2072
2073       ENDIF
2074
2075    END SUBROUTINE netcdf_handle_error_rad
[1585]2076#endif
2077
2078
[1551]2079!------------------------------------------------------------------------------!
2080! Description:
2081! ------------
[1682]2082!> Calculate temperature tendency due to radiative cooling/heating.
2083!> Cache-optimized version.
[1551]2084!------------------------------------------------------------------------------!
[1585]2085    SUBROUTINE radiation_tendency_ij ( i, j, tend )
[1496]2086
[1585]2087       USE cloud_parameters,                                                   &
[1691]2088           ONLY:  pt_d_t
[1551]2089
[1585]2090       IMPLICIT NONE
2091
[1691]2092       INTEGER(iwp) :: i, j, k !< loop indices
[1585]2093
[1691]2094       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: tend !< pt tendency term
[1585]2095
2096#if defined ( __rrtmg )
2097!
[1691]2098!--    Calculate tendency based on heating rate
[1585]2099       DO k = nzb+1, nzt+1
[1691]2100          tend(k,j,i) = tend(k,j,i) + (rad_lw_hr(k,j,i) + rad_sw_hr(k,j,i))    &
2101                                      * pt_d_t(k) * d_seconds_hour
[1585]2102       ENDDO
2103
2104#endif
2105
2106    END SUBROUTINE radiation_tendency_ij
2107
2108
[1551]2109!------------------------------------------------------------------------------!
2110! Description:
2111! ------------
[1682]2112!> Calculate temperature tendency due to radiative cooling/heating.
2113!> Vector-optimized version
[1551]2114!------------------------------------------------------------------------------!
[1585]2115    SUBROUTINE radiation_tendency ( tend )
[1551]2116
[1585]2117       USE cloud_parameters,                                                   &
[1691]2118           ONLY:  pt_d_t
[1551]2119
[1585]2120       USE indices,                                                            &
2121           ONLY:  nxl, nxr, nyn, nys
2122
2123       IMPLICIT NONE
2124
[1691]2125       INTEGER(iwp) :: i, j, k !< loop indices
[1585]2126
[1691]2127       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: tend !< pt tendency term
[1585]2128
2129#if defined ( __rrtmg )
[1691]2130!
2131!--    Calculate tendency based on heating rate
[1585]2132       DO  i = nxl, nxr
2133          DO  j = nys, nyn
2134             DO k = nzb+1, nzt+1
[1691]2135                tend(k,j,i) = tend(k,j,i) + ( rad_lw_hr(k,j,i)                 &
2136                                            +  rad_sw_hr(k,j,i) ) * pt_d_t(k)  &
2137                                            * d_seconds_hour
[1585]2138             ENDDO
2139         ENDDO
2140       ENDDO
2141#endif
2142
2143    END SUBROUTINE radiation_tendency
2144
[1496]2145 END MODULE radiation_model_mod
Note: See TracBrowser for help on using the repository browser.