source: palm/trunk/SOURCE/land_surface_model.f90 @ 1796

Last change on this file since 1796 was 1789, checked in by maronga, 9 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 85.6 KB
RevLine 
[1682]1!> @file land_surface_model.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!
[1691]16! Copyright 1997-2015 Leibniz Universitaet Hannover
[1496]17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
[1789]21!
22!
23! Former revisions:
24! -----------------
25! $Id: land_surface_model.f90 1789 2016-03-10 11:02:40Z raasch $
26!
27! 1788 2016-03-10 11:01:04Z maronga
[1788]28! Bugfix: calculate lambda_surface based on temperature gradient between skin
29! layer and soil layer instead of Obukhov length
30! Changed: moved calculation of surface specific humidity to energy balance solver
31! New: water surfaces are available by using a fixed sea surface temperature.
32! The roughness lengths are calculated dynamically using the Charnock
33! parameterization. This involves the new roughness length for moisture z0q.
34! New: modified solution of the energy balance solver and soil model for
35! paved surfaces (i.e. asphalt concrete).
36! Syntax layout improved.
37! Changed: parameter dewfall removed.
[1758]38!
[1784]39! 1783 2016-03-06 18:36:17Z raasch
40! netcdf variables moved to netcdf module
41!
[1758]42! 1757 2016-02-22 15:49:32Z maronga
[1757]43! Bugfix: set tm_soil_m to zero after allocation. Added parameter
44! unscheduled_radiation_calls to control calls of the radiation model based on
45! the skin temperature change during one time step (preliminary version). Set
46! qsws_soil_eb to zero at model start (previously set to qsws_eb). Removed MAX
47! function as it cannot be vectorized.
[1710]48!
49! 1709 2015-11-04 14:47:01Z maronga
[1709]50! Renamed pt_1 and qv_1 to pt1 and qv1.
51! Bugfix: set initial values for t_surface_p in case of restart runs
52! Bugfix: zero resistance caused crash when using radiation_scheme = 'clear-sky'
53! Bugfix: calculation of rad_net when using radiation_scheme = 'clear-sky'
54! Added todo action
[1698]55!
56! 1697 2015-10-28 17:14:10Z raasch
57! bugfix: misplaced cpp-directive
58!
59! 1695 2015-10-27 10:03:11Z maronga
[1692]60! Bugfix: REAL constants provided with KIND-attribute in call of
[1696]61! Replaced rif with ol
62!
[1692]63! 1691 2015-10-26 16:17:44Z maronga
[1691]64! Added skip_time_do_lsm to allow for spin-ups without LSM. Various bugfixes:
65! Soil temperatures are now defined at the edges of the layers, calculation of
66! shb_eb corrected, prognostic equation for skin temperature corrected. Surface
67! fluxes are now directly transfered to atmosphere
[1552]68!
[1683]69! 1682 2015-10-07 23:56:08Z knoop
70! Code annotations made doxygen readable
71!
[1591]72! 1590 2015-05-08 13:56:27Z maronga
73! Bugfix: definition of character strings requires same length for all elements
74!
[1586]75! 1585 2015-04-30 07:05:52Z maronga
76! Modifications for RRTMG. Changed tables to PARAMETER type.
77!
[1572]78! 1571 2015-03-12 16:12:49Z maronga
79! Removed upper-case variable names. Corrected distribution of precipitation to
80! the liquid water reservoir and the bare soil fractions.
81!
[1556]82! 1555 2015-03-04 17:44:27Z maronga
83! Added output of r_a and r_s
84!
[1554]85! 1553 2015-03-03 17:33:54Z maronga
86! Improved better treatment of roughness lengths. Added default soil temperature
87! profile
88!
[1552]89! 1551 2015-03-03 14:18:16Z maronga
[1551]90! Flux calculation is now done in prandtl_fluxes. Added support for data output.
91! Vertical indices have been replaced. Restart runs are now possible. Some
92! variables have beem renamed. Bugfix in the prognostic equation for the surface
93! temperature. Introduced z0_eb and z0h_eb, which overwrite the setting of
94! roughness_length and z0_factor. Added Clapp & Hornberger parametrization for
95! the hydraulic conductivity. Bugfix for root fraction and extraction
96! calculation
[1496]97!
[1514]98! intrinsic function MAX and MIN
[1496]99!
[1501]100! 1500 2014-12-03 17:42:41Z maronga
101! Corrected calculation of aerodynamic resistance (r_a).
102! Precipitation is now added to liquid water reservoir using LE_liq.
103! Added support for dry runs.
104!
[1497]105! 1496 2014-12-02 17:25:50Z maronga
106! Initial revision
107!
[1496]108!
109! Description:
110! ------------
[1682]111!> Land surface model, consisting of a solver for the energy balance at the
112!> surface and a four layer soil scheme. The scheme is similar to the TESSEL
113!> scheme implemented in the ECMWF IFS model, with modifications according to
114!> H-TESSEL. The implementation is based on the formulation implemented in the
115!> DALES and UCLA-LES models.
116!>
[1691]117!> @todo Consider partial absorption of the net shortwave radiation by the
[1709]118!>       skin layer.
[1788]119!> @todo Improve surface water parameterization
[1691]120!> @todo Invert indices (running from -3 to 0. Currently: nzb_soil=0,
121!>       nzt_soil=3)).
122!> @todo Implement surface runoff model (required when performing long-term LES
123!>       with considerable precipitation.
[1709]124!> @todo Fix crashes with radiation_scheme == 'constant'
[1682]125!>
[1709]126!> @note No time step criterion is required as long as the soil layers do not
127!>       become too thin.
[1496]128!------------------------------------------------------------------------------!
[1682]129 MODULE land_surface_model_mod
130 
[1691]131    USE arrays_3d,                                                             &
[1788]132        ONLY:  hyp, ol, pt, pt_p, q, q_p, ql, qsws, shf, ts, us, vpt, z0, z0h, &
133               z0q
[1496]134
[1691]135    USE cloud_parameters,                                                      &
136        ONLY:  cp, hyrho, l_d_cp, l_d_r, l_v, prr, pt_d_t, rho_l, r_d, r_v
[1496]137
[1691]138    USE control_parameters,                                                    &
139        ONLY:  cloud_physics, dt_3d, humidity, intermediate_timestep_count,    &
140               initializing_actions, intermediate_timestep_count_max,          &
141               max_masks, precipitation, pt_surface, rho_surface,              &
142               roughness_length, surface_pressure, timestep_scheme, tsc,       &
143               z0h_factor, time_since_reference_point
[1496]144
[1691]145    USE indices,                                                               &
146        ONLY:  nbgp, nxlg, nxrg, nyng, nysg, nzb, nzb_s_inner 
[1496]147
[1691]148    USE kinds
[1496]149
[1691]150    USE pegrid
[1496]151
[1691]152    USE radiation_model_mod,                                                   &
153        ONLY:  force_radiation_call, radiation_scheme, rad_net, rad_sw_in,     &
[1757]154               rad_lw_out, rad_lw_out_change_0, sigma_sb,                      &
155               unscheduled_radiation_calls
[1691]156       
157#if defined ( __rrtmg )
158    USE radiation_model_mod,                                                   &
[1709]159        ONLY:  rrtm_idrv
[1691]160#endif
[1496]161
[1691]162    USE statistics,                                                            &
163        ONLY:  hom, statistic_regions
164
[1496]165    IMPLICIT NONE
166
167!
168!-- LSM model constants
[1682]169    INTEGER(iwp), PARAMETER :: nzb_soil = 0, & !< bottom of the soil model (to be switched)
170                               nzt_soil = 3, & !< top of the soil model (to be switched)
171                               nzs = 4         !< number of soil layers (fixed for now)
[1496]172
173    REAL(wp), PARAMETER ::                     &
[1551]174              b_ch               = 6.04_wp,    & ! Clapp & Hornberger exponent
175              lambda_h_dry       = 0.19_wp,    & ! heat conductivity for dry soil   
[1496]176              lambda_h_sm        = 3.44_wp,    & ! heat conductivity of the soil matrix
177              lambda_h_water     = 0.57_wp,    & ! heat conductivity of water
178              psi_sat            = -0.388_wp,  & ! soil matrix potential at saturation
[1551]179              rho_c_soil         = 2.19E6_wp,  & ! volumetric heat capacity of soil
180              rho_c_water        = 4.20E6_wp,  & ! volumetric heat capacity of water
[1496]181              m_max_depth        = 0.0002_wp     ! Maximum capacity of the water reservoir (m)
182
183
184!
185!-- LSM variables
[1788]186    INTEGER(iwp) :: veg_type  = 2, & !< NAMELIST veg_type_2d
187                    soil_type = 3    !< NAMELIST soil_type_2d
[1496]188
[1788]189    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  soil_type_2d, &  !< soil type, 0: user-defined, 1-7: generic (see list)
190                                                  veg_type_2d      !< vegetation type, 0: user-defined, 1-19: generic (see list)
191
192    LOGICAL, DIMENSION(:,:), ALLOCATABLE :: water_surface,     & !< flag parameter for water surfaces (classes 14+15)
193                                            pave_surface,      & !< flag parameter for pavements (asphalt etc.) (class 20)
194                                            building_surface     !< flag parameter indicating that the surface element is covered by buildings (no LSM actions, not implemented yet)
195
[1691]196    LOGICAL :: conserve_water_content = .TRUE.,  & !< open or closed bottom surface for the soil model
197               force_radiation_call_l = .FALSE., & !< flag parameter for unscheduled radiation model calls
198               land_surface = .FALSE.              !< flag parameter indicating wheather the lsm is used
[1496]199
200!   value 9999999.9_wp -> generic available or user-defined value must be set
201!   otherwise -> no generic variable and user setting is optional
[1682]202    REAL(wp) :: alpha_vangenuchten = 9999999.9_wp,      & !< NAMELIST alpha_vg
203                canopy_resistance_coefficient = 9999999.9_wp, & !< NAMELIST g_d
204                c_surface   = 20000.0_wp,               & !< Surface (skin) heat capacity
205                drho_l_lv,                              & !< (rho_l * l_v)**-1
206                exn,                                    & !< value of the Exner function
207                e_s = 0.0_wp,                           & !< saturation water vapour pressure
208                field_capacity = 9999999.9_wp,          & !< NAMELIST m_fc
209                f_shortwave_incoming = 9999999.9_wp,    & !< NAMELIST f_sw_in
210                hydraulic_conductivity = 9999999.9_wp,  & !< NAMELIST gamma_w_sat
211                ke = 0.0_wp,                            & !< Kersten number
[1691]212                lambda_h_sat = 0.0_wp,                  & !< heat conductivity for saturated soil
[1682]213                lambda_surface_stable = 9999999.9_wp,   & !< NAMELIST lambda_surface_s
214                lambda_surface_unstable = 9999999.9_wp, & !< NAMELIST lambda_surface_u
215                leaf_area_index = 9999999.9_wp,         & !< NAMELIST lai
216                l_vangenuchten = 9999999.9_wp,          & !< NAMELIST l_vg
217                min_canopy_resistance = 9999999.9_wp,   & !< NAMELIST r_canopy_min
218                min_soil_resistance = 50.0_wp,          & !< NAMELIST r_soil_min
219                m_total = 0.0_wp,                       & !< weighted total water content of the soil (m3/m3)
220                n_vangenuchten = 9999999.9_wp,          & !< NAMELIST n_vg
[1788]221                pave_depth = 9999999.9_wp,              & !< depth of the pavement
222                pave_heat_capacity = 1.94E6_wp,         & !< volumetric heat capacity of pavement (e.g. roads)
223                pave_heat_conductivity = 1.00_wp,       & !< heat conductivity for pavements (e.g. roads)
[1682]224                q_s = 0.0_wp,                           & !< saturation specific humidity
225                residual_moisture = 9999999.9_wp,       & !< NAMELIST m_res
226                rho_cp,                                 & !< rho_surface * cp
227                rho_lv,                                 & !< rho * l_v
228                rd_d_rv,                                & !< r_d / r_v
229                saturation_moisture = 9999999.9_wp,     & !< NAMELIST m_sat
[1691]230                skip_time_do_lsm = 0.0_wp,              & !< LSM is not called before this time
[1682]231                vegetation_coverage = 9999999.9_wp,     & !< NAMELIST c_veg
232                wilting_point = 9999999.9_wp,           & !< NAMELIST m_wilt
233                z0_eb  = 9999999.9_wp,                  & !< NAMELIST z0 (lsm_par)
[1788]234                z0h_eb = 9999999.9_wp,                  & !< NAMELIST z0h (lsm_par)
235                z0q_eb = 9999999.9_wp                     !< NAMELIST z0q (lsm_par)
[1496]236
[1551]237    REAL(wp), DIMENSION(nzb_soil:nzt_soil) :: &
[1691]238              ddz_soil_stag,                  & !< 1/dz_soil_stag
239              dz_soil_stag,                   & !< soil grid spacing (center-center)
240              root_extr = 0.0_wp,             & !< root extraction
[1551]241              root_fraction = (/9999999.9_wp, 9999999.9_wp,    &
[1682]242                                9999999.9_wp, 9999999.9_wp /), & !< distribution of root surface area to the individual soil layers
243              zs = (/0.07_wp, 0.28_wp, 1.00_wp,  2.89_wp/),    & !< soil layer depths (m)
244              soil_moisture = 0.0_wp          !< soil moisture content (m3/m3)
[1496]245
[1551]246    REAL(wp), DIMENSION(nzb_soil:nzt_soil+1) ::   &
[1691]247              soil_temperature = (/290.0_wp, 287.0_wp, 285.0_wp,  283.0_wp,    & !< soil temperature (K)
248                                   283.0_wp /),                                &                                   
249              ddz_soil,                                                        & !< 1/dz_soil
250              dz_soil                                                            !< soil grid spacing (edge-edge)
[1496]251
252#if defined( __nopointer )
[1682]253    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET :: t_surface,   & !< surface temperature (K)
254                                                     t_surface_p, & !< progn. surface temperature (K)
255                                                     m_liq_eb,    & !< liquid water reservoir (m)
256                                                     m_liq_eb_av, & !< liquid water reservoir (m)
257                                                     m_liq_eb_p     !< progn. liquid water reservoir (m)
[1496]258#else
[1551]259    REAL(wp), DIMENSION(:,:), POINTER :: t_surface,      &
260                                         t_surface_p,    & 
261                                         m_liq_eb,       & 
262                                         m_liq_eb_p
[1496]263
[1551]264    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET :: t_surface_1, t_surface_2, &
265                                                     m_liq_eb_av,              &
266                                                     m_liq_eb_1, m_liq_eb_2
[1496]267#endif
268
269!
270!-- Temporal tendencies for time stepping           
[1682]271    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: tt_surface_m,  & !< surface temperature tendency (K)
272                                             tm_liq_eb_m      !< liquid water reservoir tendency (m)
[1496]273
274!
275!-- Energy balance variables               
[1691]276    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: &
[1682]277              alpha_vg,         & !< coef. of Van Genuchten
278              c_liq,            & !< liquid water coverage (of vegetated area)
279              c_liq_av,         & !< average of c_liq
280              c_soil_av,        & !< average of c_soil
281              c_veg,            & !< vegetation coverage
282              c_veg_av,         & !< average of c_veg
283              f_sw_in,          & !< fraction of absorbed shortwave radiation by the surface layer (not implemented yet)
284              ghf_eb,           & !< ground heat flux
285              ghf_eb_av,        & !< average of ghf_eb
286              gamma_w_sat,      & !< hydraulic conductivity at saturation
287              g_d,              & !< coefficient for dependence of r_canopy on water vapour pressure deficit
288              lai,              & !< leaf area index
289              lai_av,           & !< average of lai
[1691]290              lambda_surface_s, & !< coupling between surface and soil (depends on vegetation type)
291              lambda_surface_u, & !< coupling between surface and soil (depends on vegetation type)
[1682]292              l_vg,             & !< coef. of Van Genuchten
293              m_fc,             & !< soil moisture at field capacity (m3/m3)
294              m_res,            & !< residual soil moisture
295              m_sat,            & !< saturation soil moisture (m3/m3)
296              m_wilt,           & !< soil moisture at permanent wilting point (m3/m3)
297              n_vg,             & !< coef. Van Genuchten 
298              qsws_eb,          & !< surface flux of latent heat (total)
299              qsws_eb_av,       & !< average of qsws_eb
300              qsws_liq_eb,      & !< surface flux of latent heat (liquid water portion)
301              qsws_liq_eb_av,   & !< average of qsws_liq_eb
302              qsws_soil_eb,     & !< surface flux of latent heat (soil portion)
303              qsws_soil_eb_av,  & !< average of qsws_soil_eb
304              qsws_veg_eb,      & !< surface flux of latent heat (vegetation portion)
305              qsws_veg_eb_av,   & !< average of qsws_veg_eb
306              rad_net_l,        & !< local copy of rad_net (net radiation at surface)
307              r_a,              & !< aerodynamic resistance
[1709]308              r_a_av,           & !< average of r_a
[1682]309              r_canopy,         & !< canopy resistance
[1691]310              r_soil,           & !< soil resistance
[1682]311              r_soil_min,       & !< minimum soil resistance
312              r_s,              & !< total surface resistance (combination of r_soil and r_canopy)
[1691]313              r_s_av,           & !< average of r_s
[1682]314              r_canopy_min,     & !< minimum canopy (stomatal) resistance
315              shf_eb,           & !< surface flux of sensible heat
316              shf_eb_av           !< average of shf_eb
[1496]317
[1691]318
[1496]319    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::                                 &
[1788]320              lambda_h, &   !< heat conductivity of soil (W/m/K)                           
[1682]321              lambda_w, &   !< hydraulic diffusivity of soil (?)
[1788]322              gamma_w,  &   !< hydraulic conductivity of soil (W/m/K)
[1682]323              rho_c_total   !< volumetric heat capacity of the actual soil matrix (?)
[1496]324
325#if defined( __nopointer )
326    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                         &
[1682]327              t_soil,    & !< Soil temperature (K)
328              t_soil_av, & !< Average of t_soil
329              t_soil_p,  & !< Prog. soil temperature (K)
330              m_soil,    & !< Soil moisture (m3/m3)
331              m_soil_av, & !< Average of m_soil
332              m_soil_p     !< Prog. soil moisture (m3/m3)
[1496]333#else
334    REAL(wp), DIMENSION(:,:,:), POINTER ::                                     &
[1551]335              t_soil, t_soil_p, &
[1496]336              m_soil, m_soil_p   
337
338    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                         &
[1551]339              t_soil_av, t_soil_1, t_soil_2,                                   &
340              m_soil_av, m_soil_1, m_soil_2
[1496]341#endif
342
343
344    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::                                 &
[1757]345              tt_soil_m, & !< t_soil storage array
346              tm_soil_m, & !< m_soil storage array
[1682]347              root_fr      !< root fraction (sum=1)
[1496]348
[1551]349
[1496]350!
[1551]351!-- Predefined Land surface classes (veg_type)
[1788]352    CHARACTER(26), DIMENSION(0:20), PARAMETER :: veg_type_name = (/ &
[1590]353                                   'user defined              ',    & ! 0
354                                   'crops, mixed farming      ',    & !  1
355                                   'short grass               ',    & !  2
[1585]356                                   'evergreen needleleaf trees',    & !  3
357                                   'deciduous needleleaf trees',    & !  4
[1590]358                                   'evergreen broadleaf trees ',    & !  5
359                                   'deciduous broadleaf trees ',    & !  6
360                                   'tall grass                ',    & !  7
361                                   'desert                    ',    & !  8
362                                   'tundra                    ',    & !  9
363                                   'irrigated crops           ',    & ! 10
364                                   'semidesert                ',    & ! 11
365                                   'ice caps and glaciers     ',    & ! 12
366                                   'bogs and marshes          ',    & ! 13
367                                   'inland water              ',    & ! 14
368                                   'ocean                     ',    & ! 15
369                                   'evergreen shrubs          ',    & ! 16
370                                   'deciduous shrubs          ',    & ! 17
371                                   'mixed forest/woodland     ',    & ! 18
[1788]372                                   'interrupted forest        ',    & ! 19
373                                   'pavements/roads           '     & ! 20
[1585]374                                                                 /)
[1496]375
376!
[1551]377!-- Soil model classes (soil_type)
[1585]378    CHARACTER(12), DIMENSION(0:7), PARAMETER :: soil_type_name = (/ &
379                                   'user defined',                  & ! 0
[1590]380                                   'coarse      ',                  & ! 1
381                                   'medium      ',                  & ! 2
382                                   'medium-fine ',                  & ! 3
383                                   'fine        ',                  & ! 4
384                                   'very fine   ',                  & ! 5
385                                   'organic     ',                  & ! 6
386                                   'loamy (CH)  '                   & ! 7
[1585]387                                                                 /)
[1551]388!
389!-- Land surface parameters according to the respective classes (veg_type)
390
391!
392!-- Land surface parameters I
393!--                          r_canopy_min,     lai,   c_veg,     g_d
[1788]394    REAL(wp), DIMENSION(0:3,1:20), PARAMETER :: veg_pars = RESHAPE( (/ &
[1585]395                                 180.0_wp, 3.00_wp, 0.90_wp, 0.00_wp,  & !  1
396                                 110.0_wp, 2.00_wp, 0.85_wp, 0.00_wp,  & !  2
397                                 500.0_wp, 5.00_wp, 0.90_wp, 0.03_wp,  & !  3
398                                 500.0_wp, 5.00_wp, 0.90_wp, 0.03_wp,  & !  4
399                                 175.0_wp, 5.00_wp, 0.90_wp, 0.03_wp,  & !  5
400                                 240.0_wp, 6.00_wp, 0.99_wp, 0.13_wp,  & !  6
401                                 100.0_wp, 2.00_wp, 0.70_wp, 0.00_wp,  & !  7
402                                 250.0_wp, 0.05_wp, 0.00_wp, 0.00_wp,  & !  8
403                                  80.0_wp, 1.00_wp, 0.50_wp, 0.00_wp,  & !  9
404                                 180.0_wp, 3.00_wp, 0.90_wp, 0.00_wp,  & ! 10
405                                 150.0_wp, 0.50_wp, 0.10_wp, 0.00_wp,  & ! 11
406                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp,  & ! 12
407                                 240.0_wp, 4.00_wp, 0.60_wp, 0.00_wp,  & ! 13
408                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp,  & ! 14
409                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp,  & ! 15
410                                 225.0_wp, 3.00_wp, 0.50_wp, 0.00_wp,  & ! 16
411                                 225.0_wp, 1.50_wp, 0.50_wp, 0.00_wp,  & ! 17
412                                 250.0_wp, 5.00_wp, 0.90_wp, 0.03_wp,  & ! 18
[1788]413                                 175.0_wp, 2.50_wp, 0.90_wp, 0.03_wp,  & ! 19
414                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp   & ! 20
415                                 /), (/ 4, 20 /) )
[1496]416
417!
[1788]418!-- Land surface parameters II          z0,         z0h,         z0q
419    REAL(wp), DIMENSION(0:2,1:20), PARAMETER :: roughness_par = RESHAPE( (/ & 
420                                   0.25_wp,  0.25E-2_wp,  0.25E-2_wp,       & !  1
421                                   0.20_wp,  0.20E-2_wp,  0.20E-2_wp,       & !  2
422                                   2.00_wp,     2.00_wp,     2.00_wp,       & !  3
423                                   2.00_wp,     2.00_wp,     2.00_wp,       & !  4
424                                   2.00_wp,     2.00_wp,     2.00_wp,       & !  5
425                                   2.00_wp,     2.00_wp,     2.00_wp,       & !  6
426                                   0.47_wp,  0.47E-2_wp,  0.47E-2_wp,       & !  7
427                                  0.013_wp, 0.013E-2_wp, 0.013E-2_wp,       & !  8
428                                  0.034_wp, 0.034E-2_wp, 0.034E-2_wp,       & !  9
429                                    0.5_wp,  0.50E-2_wp,  0.50E-2_wp,       & ! 10
430                                   0.17_wp,  0.17E-2_wp,  0.17E-2_wp,       & ! 11
431                                 1.3E-3_wp,   1.3E-4_wp,   1.3E-4_wp,       & ! 12
432                                   0.83_wp,  0.83E-2_wp,  0.83E-2_wp,       & ! 13
433                                   0.00_wp,     0.00_wp,     0.00_wp,       & ! 14
434                                   0.00_wp,     0.00_wp,     0.00_wp,       & ! 15
435                                   0.10_wp,  0.10E-2_wp,  0.10E-2_wp,       & ! 16
436                                   0.25_wp,  0.25E-2_wp,  0.25E-2_wp,       & ! 17
437                                   2.00_wp,  2.00E-2_wp,  2.00E-2_wp,       & ! 18
438                                   1.10_wp,  1.10E-2_wp,  1.10E-2_wp,       & ! 19
439                                 1.0E-4_wp,   1.0E-5_wp,   1.0E-5_wp        & ! 20
440                                 /), (/ 3, 20 /) )
[1496]441
442!
[1551]443!-- Land surface parameters III lambda_surface_s, lambda_surface_u, f_sw_in
[1788]444    REAL(wp), DIMENSION(0:2,1:20), PARAMETER :: surface_pars = RESHAPE( (/ &
[1585]445                                      10.0_wp,       10.0_wp, 0.05_wp,     & !  1
446                                      10.0_wp,       10.0_wp, 0.05_wp,     & !  2
447                                      20.0_wp,       15.0_wp, 0.03_wp,     & !  3
448                                      20.0_wp,       15.0_wp, 0.03_wp,     & !  4
449                                      20.0_wp,       15.0_wp, 0.03_wp,     & !  5
450                                      20.0_wp,       15.0_wp, 0.03_wp,     & !  6
451                                      10.0_wp,       10.0_wp, 0.05_wp,     & !  7
452                                      15.0_wp,       15.0_wp, 0.00_wp,     & !  8
453                                      10.0_wp,       10.0_wp, 0.05_wp,     & !  9
454                                      10.0_wp,       10.0_wp, 0.05_wp,     & ! 10
455                                      10.0_wp,       10.0_wp, 0.05_wp,     & ! 11
456                                      58.0_wp,       58.0_wp, 0.00_wp,     & ! 12
457                                      10.0_wp,       10.0_wp, 0.05_wp,     & ! 13
[1788]458                                    1.0E10_wp,     1.0E10_wp, 0.00_wp,     & ! 14
459                                    1.0E10_wp,     1.0E10_wp, 0.00_wp,     & ! 15
[1585]460                                      10.0_wp,       10.0_wp, 0.05_wp,     & ! 16
461                                      10.0_wp,       10.0_wp, 0.05_wp,     & ! 17
462                                      20.0_wp,       15.0_wp, 0.03_wp,     & ! 18
[1788]463                                      20.0_wp,       15.0_wp, 0.03_wp,     & ! 19
464                                       0.0_wp,        0.0_wp, 0.00_wp      & ! 20
465                                      /), (/ 3, 20 /) )
[1496]466
467!
468!-- Root distribution (sum = 1)  level 1, level 2, level 3, level 4,
[1788]469    REAL(wp), DIMENSION(0:3,1:20), PARAMETER :: root_distribution = RESHAPE( (/ &
[1585]470                                 0.24_wp, 0.41_wp, 0.31_wp, 0.04_wp,            & !  1
471                                 0.35_wp, 0.38_wp, 0.23_wp, 0.04_wp,            & !  2
472                                 0.26_wp, 0.39_wp, 0.29_wp, 0.06_wp,            & !  3
473                                 0.26_wp, 0.38_wp, 0.29_wp, 0.07_wp,            & !  4
474                                 0.24_wp, 0.38_wp, 0.31_wp, 0.07_wp,            & !  5
475                                 0.25_wp, 0.34_wp, 0.27_wp, 0.14_wp,            & !  6
476                                 0.27_wp, 0.27_wp, 0.27_wp, 0.09_wp,            & !  7
477                                 1.00_wp, 0.00_wp, 0.00_wp, 0.00_wp,            & !  8
478                                 0.47_wp, 0.45_wp, 0.08_wp, 0.00_wp,            & !  9
479                                 0.24_wp, 0.41_wp, 0.31_wp, 0.04_wp,            & ! 10
480                                 0.17_wp, 0.31_wp, 0.33_wp, 0.19_wp,            & ! 11
481                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp,            & ! 12
482                                 0.25_wp, 0.34_wp, 0.27_wp, 0.11_wp,            & ! 13
483                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp,            & ! 14
484                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp,            & ! 15
485                                 0.23_wp, 0.36_wp, 0.30_wp, 0.11_wp,            & ! 16
486                                 0.23_wp, 0.36_wp, 0.30_wp, 0.11_wp,            & ! 17
487                                 0.19_wp, 0.35_wp, 0.36_wp, 0.10_wp,            & ! 18
[1788]488                                 0.19_wp, 0.35_wp, 0.36_wp, 0.10_wp,            & ! 19
489                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp             & ! 20
490                                 /), (/ 4, 20 /) )
[1496]491
492!
493!-- Soil parameters according to the following porosity classes (soil_type)
[1551]494
[1496]495!
[1551]496!-- Soil parameters I           alpha_vg,      l_vg,    n_vg, gamma_w_sat
[1585]497    REAL(wp), DIMENSION(0:3,1:7), PARAMETER :: soil_pars = RESHAPE( (/     &
[1496]498                                 3.83_wp,  1.250_wp, 1.38_wp,  6.94E-6_wp, & ! 1
499                                 3.14_wp, -2.342_wp, 1.28_wp,  1.16E-6_wp, & ! 2
500                                 0.83_wp, -0.588_wp, 1.25_wp,  0.26E-6_wp, & ! 3
501                                 3.67_wp, -1.977_wp, 1.10_wp,  2.87E-6_wp, & ! 4
502                                 2.65_wp,  2.500_wp, 1.10_wp,  1.74E-6_wp, & ! 5
[1551]503                                 1.30_wp,  0.400_wp, 1.20_wp,  0.93E-6_wp, & ! 6
504                                 0.00_wp,  0.00_wp,  0.00_wp,  0.57E-6_wp  & ! 7
505                                 /), (/ 4, 7 /) )
[1496]506
507!
508!-- Soil parameters II              m_sat,     m_fc,   m_wilt,    m_res 
[1585]509    REAL(wp), DIMENSION(0:3,1:7), PARAMETER :: m_soil_pars = RESHAPE( (/ &
[1496]510                                 0.403_wp, 0.244_wp, 0.059_wp, 0.025_wp, & ! 1
511                                 0.439_wp, 0.347_wp, 0.151_wp, 0.010_wp, & ! 2
512                                 0.430_wp, 0.383_wp, 0.133_wp, 0.010_wp, & ! 3
513                                 0.520_wp, 0.448_wp, 0.279_wp, 0.010_wp, & ! 4
514                                 0.614_wp, 0.541_wp, 0.335_wp, 0.010_wp, & ! 5
[1551]515                                 0.766_wp, 0.663_wp, 0.267_wp, 0.010_wp, & ! 6
516                                 0.472_wp, 0.323_wp, 0.171_wp, 0.000_wp  & ! 7
517                                 /), (/ 4, 7 /) )
[1496]518
519
520    SAVE
521
522
523    PRIVATE
524
525
[1551]526!
527!-- Public parameters, constants and initial values
528    PUBLIC alpha_vangenuchten, c_surface, canopy_resistance_coefficient,       &
[1788]529           conserve_water_content, field_capacity,                             &
[1551]530           f_shortwave_incoming, hydraulic_conductivity, init_lsm,             &
531           init_lsm_arrays, lambda_surface_stable, lambda_surface_unstable,    &
532           land_surface, leaf_area_index, lsm_energy_balance, lsm_soil_model,  &
533           lsm_swap_timelevel, l_vangenuchten, min_canopy_resistance,          &
[1788]534           min_soil_resistance, n_vangenuchten, pave_heat_capacity,            &
535           pave_depth, pave_heat_conductivity, residual_moisture, rho_cp,      &
[1691]536           rho_lv, root_fraction, saturation_moisture, skip_time_do_lsm,       &
537           soil_moisture, soil_temperature, soil_type, soil_type_name,         &
538           vegetation_coverage, veg_type, veg_type_name, wilting_point, z0_eb, &
[1788]539           z0h_eb, z0q_eb
[1496]540
[1551]541!
[1783]542!-- Public grid variables
543    PUBLIC nzb_soil, nzs, nzt_soil, zs
[1496]544
[1551]545!
546!-- Public 2D output variables
547    PUBLIC c_liq, c_liq_av, c_soil_av, c_veg, c_veg_av, ghf_eb, ghf_eb_av,     &
548           lai, lai_av, qsws_eb, qsws_eb_av, qsws_liq_eb, qsws_liq_eb_av,      &
549           qsws_soil_eb, qsws_soil_eb_av, qsws_veg_eb, qsws_veg_eb_av,         &
[1555]550           r_a, r_a_av, r_s, r_s_av, shf_eb, shf_eb_av
[1551]551
552!
553!-- Public prognostic variables
554    PUBLIC m_liq_eb, m_liq_eb_av, m_soil, m_soil_av, t_soil, t_soil_av
555
[1496]556    INTERFACE init_lsm
557       MODULE PROCEDURE init_lsm
558    END INTERFACE init_lsm
559
560    INTERFACE lsm_energy_balance
561       MODULE PROCEDURE lsm_energy_balance
562    END INTERFACE lsm_energy_balance
563
564    INTERFACE lsm_soil_model
565       MODULE PROCEDURE lsm_soil_model
566    END INTERFACE lsm_soil_model
567
[1551]568    INTERFACE lsm_swap_timelevel
569       MODULE PROCEDURE lsm_swap_timelevel
570    END INTERFACE lsm_swap_timelevel
[1496]571
572 CONTAINS
573
574
575!------------------------------------------------------------------------------!
576! Description:
577! ------------
[1682]578!> Allocate land surface model arrays and define pointers
[1496]579!------------------------------------------------------------------------------!
[1551]580    SUBROUTINE init_lsm_arrays
[1496]581   
582
583       IMPLICIT NONE
584
585!
[1551]586!--    Allocate surface and soil temperature / humidity
[1496]587#if defined( __nopointer )
[1551]588       ALLOCATE ( m_liq_eb(nysg:nyng,nxlg:nxrg) )
589       ALLOCATE ( m_liq_eb_p(nysg:nyng,nxlg:nxrg) )
590       ALLOCATE ( m_soil(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
591       ALLOCATE ( m_soil_p(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
592       ALLOCATE ( t_surface(nysg:nyng,nxlg:nxrg) )
593       ALLOCATE ( t_surface_p(nysg:nyng,nxlg:nxrg) )
594       ALLOCATE ( t_soil(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
595       ALLOCATE ( t_soil_p(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
[1496]596#else
[1551]597       ALLOCATE ( m_liq_eb_1(nysg:nyng,nxlg:nxrg) )
598       ALLOCATE ( m_liq_eb_2(nysg:nyng,nxlg:nxrg) )
599       ALLOCATE ( m_soil_1(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
600       ALLOCATE ( m_soil_2(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
601       ALLOCATE ( t_surface_1(nysg:nyng,nxlg:nxrg) )
602       ALLOCATE ( t_surface_2(nysg:nyng,nxlg:nxrg) )
603       ALLOCATE ( t_soil_1(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
604       ALLOCATE ( t_soil_2(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
[1496]605#endif
606
607!
[1551]608!--    Allocate intermediate timestep arrays
609       ALLOCATE ( tm_liq_eb_m(nysg:nyng,nxlg:nxrg) )
610       ALLOCATE ( tm_soil_m(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
611       ALLOCATE ( tt_surface_m(nysg:nyng,nxlg:nxrg) )
612       ALLOCATE ( tt_soil_m(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
[1496]613
614!
615!--    Allocate 2D vegetation model arrays
[1551]616       ALLOCATE ( alpha_vg(nysg:nyng,nxlg:nxrg) )
[1788]617       ALLOCATE ( building_surface(nysg:nyng,nxlg:nxrg) )
[1496]618       ALLOCATE ( c_liq(nysg:nyng,nxlg:nxrg) )
619       ALLOCATE ( c_veg(nysg:nyng,nxlg:nxrg) )
[1551]620       ALLOCATE ( f_sw_in(nysg:nyng,nxlg:nxrg) )
621       ALLOCATE ( ghf_eb(nysg:nyng,nxlg:nxrg) )
[1496]622       ALLOCATE ( gamma_w_sat(nysg:nyng,nxlg:nxrg) )
[1551]623       ALLOCATE ( g_d(nysg:nyng,nxlg:nxrg) )
624       ALLOCATE ( lai(nysg:nyng,nxlg:nxrg) )
625       ALLOCATE ( l_vg(nysg:nyng,nxlg:nxrg) )
626       ALLOCATE ( lambda_surface_u(nysg:nyng,nxlg:nxrg) )
627       ALLOCATE ( lambda_surface_s(nysg:nyng,nxlg:nxrg) )
[1496]628       ALLOCATE ( m_fc(nysg:nyng,nxlg:nxrg) )
629       ALLOCATE ( m_res(nysg:nyng,nxlg:nxrg) )
630       ALLOCATE ( m_sat(nysg:nyng,nxlg:nxrg) )
631       ALLOCATE ( m_wilt(nysg:nyng,nxlg:nxrg) )
[1551]632       ALLOCATE ( n_vg(nysg:nyng,nxlg:nxrg) )
[1788]633       ALLOCATE ( pave_surface(nysg:nyng,nxlg:nxrg) )
[1551]634       ALLOCATE ( qsws_eb(nysg:nyng,nxlg:nxrg) )
635       ALLOCATE ( qsws_soil_eb(nysg:nyng,nxlg:nxrg) )
636       ALLOCATE ( qsws_liq_eb(nysg:nyng,nxlg:nxrg) )
637       ALLOCATE ( qsws_veg_eb(nysg:nyng,nxlg:nxrg) )
[1585]638       ALLOCATE ( rad_net_l(nysg:nyng,nxlg:nxrg) )
[1496]639       ALLOCATE ( r_a(nysg:nyng,nxlg:nxrg) )
640       ALLOCATE ( r_canopy(nysg:nyng,nxlg:nxrg) )
641       ALLOCATE ( r_soil(nysg:nyng,nxlg:nxrg) )
642       ALLOCATE ( r_soil_min(nysg:nyng,nxlg:nxrg) )
643       ALLOCATE ( r_s(nysg:nyng,nxlg:nxrg) )
[1551]644       ALLOCATE ( r_canopy_min(nysg:nyng,nxlg:nxrg) )
645       ALLOCATE ( shf_eb(nysg:nyng,nxlg:nxrg) )
[1788]646       ALLOCATE ( soil_type_2d(nysg:nyng,nxlg:nxrg) )
647       ALLOCATE ( veg_type_2d(nysg:nyng,nxlg:nxrg) )
648       ALLOCATE ( water_surface(nysg:nyng,nxlg:nxrg) )
[1496]649
[1551]650#if ! defined( __nopointer )
[1496]651!
[1585]652!--    Initial assignment of the pointers
[1551]653       t_soil    => t_soil_1;    t_soil_p    => t_soil_2
654       t_surface => t_surface_1; t_surface_p => t_surface_2
655       m_soil    => m_soil_1;    m_soil_p    => m_soil_2
656       m_liq_eb  => m_liq_eb_1;  m_liq_eb_p  => m_liq_eb_2
657#endif
[1496]658
659
[1551]660    END SUBROUTINE init_lsm_arrays
[1500]661
[1551]662!------------------------------------------------------------------------------!
663! Description:
664! ------------
[1682]665!> Initialization of the land surface model
[1551]666!------------------------------------------------------------------------------!
667    SUBROUTINE init_lsm
668   
669
670       IMPLICIT NONE
671
[1682]672       INTEGER(iwp) ::  i !< running index
673       INTEGER(iwp) ::  j !< running index
674       INTEGER(iwp) ::  k !< running index
[1551]675
[1709]676       REAL(wp) :: pt1   !< potential temperature at first grid level
[1551]677
[1691]678
[1551]679!
680!--    Calculate Exner function
[1691]681       exn = ( surface_pressure / 1000.0_wp )**0.286_wp
[1551]682
683
684!
685!--    If no cloud physics is used, rho_surface has not been calculated before
[1788]686       IF (  .NOT.  cloud_physics )  THEN
[1551]687          rho_surface = surface_pressure * 100.0_wp / ( r_d * pt_surface * exn )
688       ENDIF
689
690!
691!--    Calculate frequently used parameters
692       rho_cp    = cp * rho_surface
693       rd_d_rv   = r_d / r_v
694       rho_lv    = rho_surface * l_v
695       drho_l_lv = 1.0_wp / (rho_l * l_v)
696
697!
698!--    Set inital values for prognostic quantities
699       tt_surface_m = 0.0_wp
700       tt_soil_m    = 0.0_wp
[1757]701       tm_soil_m    = 0.0_wp
[1551]702       tm_liq_eb_m  = 0.0_wp
703       c_liq        = 0.0_wp
704
705       ghf_eb = 0.0_wp
706       shf_eb = rho_cp * shf
707
[1500]708       IF ( humidity )  THEN
[1551]709          qsws_eb = rho_l * l_v * qsws
[1500]710       ELSE
[1551]711          qsws_eb = 0.0_wp
[1500]712       ENDIF
713
[1551]714       qsws_liq_eb  = 0.0_wp
[1757]715       qsws_soil_eb = 0.0_wp
[1551]716       qsws_veg_eb  = 0.0_wp
[1496]717
718       r_a        = 50.0_wp
[1551]719       r_s        = 50.0_wp
[1496]720       r_canopy   = 0.0_wp
721       r_soil     = 0.0_wp
722
723!
724!--    Allocate 3D soil model arrays
[1551]725       ALLOCATE ( root_fr(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
726       ALLOCATE ( lambda_h(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
727       ALLOCATE ( rho_c_total(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
[1496]728
729       lambda_h = 0.0_wp
730!
731!--    If required, allocate humidity-related variables for the soil model
732       IF ( humidity )  THEN
[1551]733          ALLOCATE ( lambda_w(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
734          ALLOCATE ( gamma_w(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )   
[1496]735
736          lambda_w = 0.0_wp 
737       ENDIF
738
739!
740!--    Calculate grid spacings. Temperature and moisture are defined at
[1691]741!--    the edges of the soil layers (_stag), whereas gradients/fluxes are defined
742!--    at the centers
743       dz_soil(nzb_soil) = zs(nzb_soil)
[1496]744
[1691]745       DO  k = nzb_soil+1, nzt_soil
746          dz_soil(k) = zs(k) - zs(k-1)
[1496]747       ENDDO
[1691]748       dz_soil(nzt_soil+1) = dz_soil(nzt_soil)
[1496]749
[1691]750       DO  k = nzb_soil, nzt_soil-1
751          dz_soil_stag(k) = 0.5_wp * (dz_soil(k+1) + dz_soil(k))
[1496]752       ENDDO
[1691]753       dz_soil_stag(nzt_soil) = dz_soil(nzt_soil)
[1496]754
[1551]755       ddz_soil      = 1.0_wp / dz_soil
756       ddz_soil_stag = 1.0_wp / dz_soil_stag
757
[1496]758!
[1551]759!--    Initialize standard soil types. It is possible to overwrite each
760!--    parameter by setting the respecticy NAMELIST variable to a
761!--    value /= 9999999.9.
762       IF ( soil_type /= 0 )  THEN 
763 
764          IF ( alpha_vangenuchten == 9999999.9_wp )  THEN
765             alpha_vangenuchten = soil_pars(0,soil_type)
766          ENDIF
767
768          IF ( l_vangenuchten == 9999999.9_wp )  THEN
769             l_vangenuchten = soil_pars(1,soil_type)
770          ENDIF
771
772          IF ( n_vangenuchten == 9999999.9_wp )  THEN
773             n_vangenuchten = soil_pars(2,soil_type)           
774          ENDIF
775
776          IF ( hydraulic_conductivity == 9999999.9_wp )  THEN
777             hydraulic_conductivity = soil_pars(3,soil_type)           
778          ENDIF
779
780          IF ( saturation_moisture == 9999999.9_wp )  THEN
781             saturation_moisture = m_soil_pars(0,soil_type)           
782          ENDIF
783
784          IF ( field_capacity == 9999999.9_wp )  THEN
785             field_capacity = m_soil_pars(1,soil_type)           
786          ENDIF
787
788          IF ( wilting_point == 9999999.9_wp )  THEN
789             wilting_point = m_soil_pars(2,soil_type)           
790          ENDIF
791
792          IF ( residual_moisture == 9999999.9_wp )  THEN
793             residual_moisture = m_soil_pars(3,soil_type)       
794          ENDIF
795
[1496]796       ENDIF   
797
[1788]798!
799!--    Map values to the respective 2D arrays
[1551]800       alpha_vg      = alpha_vangenuchten
801       l_vg          = l_vangenuchten
802       n_vg          = n_vangenuchten 
803       gamma_w_sat   = hydraulic_conductivity
804       m_sat         = saturation_moisture
805       m_fc          = field_capacity
806       m_wilt        = wilting_point
807       m_res         = residual_moisture
808       r_soil_min    = min_soil_resistance
809
[1496]810!
[1551]811!--    Initial run actions
812       IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
[1496]813
[1551]814          t_soil    = 0.0_wp
815          m_liq_eb  = 0.0_wp
816          m_soil    = 0.0_wp
[1496]817
[1551]818!
819!--       Map user settings of T and q for each soil layer
820!--       (make sure that the soil moisture does not drop below the permanent
821!--       wilting point) -> problems with devision by zero)
[1691]822          DO  k = nzb_soil, nzt_soil
[1551]823             t_soil(k,:,:)    = soil_temperature(k)
824             m_soil(k,:,:)    = MAX(soil_moisture(k),m_wilt(:,:))
825             soil_moisture(k) = MAX(soil_moisture(k),wilting_point)
826          ENDDO
827          t_soil(nzt_soil+1,:,:) = soil_temperature(nzt_soil+1)
[1496]828
[1551]829!
830!--       Calculate surface temperature
[1691]831          t_surface   = pt_surface * exn
[1496]832
833!
[1788]834!--       Set artifical values for ts and us so that r_a has its initial value
835!--       for the first time step
[1551]836          DO  i = nxlg, nxrg
837             DO  j = nysg, nyng
838                k = nzb_s_inner(j,i)
839
[1691]840                IF ( cloud_physics )  THEN
[1709]841                   pt1 = pt(k+1,j,i) + l_d_cp * pt_d_t(k+1) * ql(k+1,j,i)
[1691]842                ELSE
[1709]843                   pt1 = pt(k+1,j,i)
[1691]844                ENDIF
845
[1551]846!
847!--             Assure that r_a cannot be zero at model start
[1709]848                IF ( pt1 == pt(k,j,i) )  pt1 = pt1 + 1.0E-10_wp
[1551]849
[1691]850                us(j,i)  = 0.1_wp
[1709]851                ts(j,i)  = (pt1 - pt(k,j,i)) / r_a(j,i)
[1551]852                shf(j,i) = - us(j,i) * ts(j,i)
853             ENDDO
854          ENDDO
855
856!
857!--    Actions for restart runs
858       ELSE
859
860          DO  i = nxlg, nxrg
861             DO  j = nysg, nyng
862                k = nzb_s_inner(j,i)               
863                t_surface(j,i) = pt(k,j,i) * exn
864             ENDDO
865          ENDDO
866
867       ENDIF
868
[1691]869       DO  k = nzb_soil, nzt_soil
[1551]870          root_fr(k,:,:) = root_fraction(k)
871       ENDDO
[1496]872
[1551]873       IF ( veg_type /= 0 )  THEN
874          IF ( min_canopy_resistance == 9999999.9_wp )  THEN
875             min_canopy_resistance = veg_pars(0,veg_type)
876          ENDIF
877          IF ( leaf_area_index == 9999999.9_wp )  THEN
878             leaf_area_index = veg_pars(1,veg_type)         
879          ENDIF
880          IF ( vegetation_coverage == 9999999.9_wp )  THEN
881             vegetation_coverage = veg_pars(2,veg_type)     
882          ENDIF
883          IF ( canopy_resistance_coefficient == 9999999.9_wp )  THEN
884              canopy_resistance_coefficient= veg_pars(3,veg_type)     
885          ENDIF
886          IF ( lambda_surface_stable == 9999999.9_wp )  THEN
887             lambda_surface_stable = surface_pars(0,veg_type)         
888          ENDIF
889          IF ( lambda_surface_unstable == 9999999.9_wp )  THEN
890             lambda_surface_unstable = surface_pars(1,veg_type)       
891          ENDIF
892          IF ( f_shortwave_incoming == 9999999.9_wp )  THEN
893             f_shortwave_incoming = surface_pars(2,veg_type)       
894          ENDIF
895          IF ( z0_eb == 9999999.9_wp )  THEN
896             roughness_length = roughness_par(0,veg_type) 
897             z0_eb            = roughness_par(0,veg_type) 
898          ENDIF
899          IF ( z0h_eb == 9999999.9_wp )  THEN
900             z0h_eb = roughness_par(1,veg_type)
901          ENDIF
[1788]902          IF ( z0q_eb == 9999999.9_wp )  THEN
903             z0q_eb = roughness_par(2,veg_type)
904          ENDIF
905          z0h_factor = z0h_eb / ( z0_eb + 1.0E-20_wp )
[1496]906
[1551]907          IF ( ANY( root_fraction == 9999999.9_wp ) )  THEN
[1691]908             DO  k = nzb_soil, nzt_soil
[1551]909                root_fr(k,:,:) = root_distribution(k,veg_type)
910                root_fraction(k) = root_distribution(k,veg_type)
911             ENDDO
912          ENDIF
[1496]913
[1553]914       ELSE
915
916          IF ( z0_eb == 9999999.9_wp )  THEN
917             z0_eb = roughness_length
918          ENDIF
919          IF ( z0h_eb == 9999999.9_wp )  THEN
920             z0h_eb = z0_eb * z0h_factor
921          ENDIF
[1788]922          IF ( z0q_eb == 9999999.9_wp )  THEN
923             z0q_eb = z0_eb * z0h_factor
924          ENDIF
[1553]925
[1496]926       ENDIF
927
928!
[1788]929!--    For surfaces covered with pavement, set depth of the pavement (with dry
930!--    soil below). The depth must be greater than the first soil layer depth
931       IF ( veg_type == 20 )  THEN
932          IF ( pave_depth == 9999999.9_wp )  THEN
933             pave_depth = zs(nzb_soil) 
934          ELSE
935             pave_depth = MAX( zs(nzb_soil), pave_depth )
936          ENDIF
937       ENDIF
938
939!
940!--    Map vegetation and soil types to 2D array to allow for heterogeneous
941!--    surfaces via user interface see below
942       veg_type_2d = veg_type
943       soil_type_2d = soil_type
944
945!
946!--    Map vegetation parameters to the respective 2D arrays
[1551]947       r_canopy_min         = min_canopy_resistance
948       lai                  = leaf_area_index
949       c_veg                = vegetation_coverage
950       g_d                  = canopy_resistance_coefficient
951       lambda_surface_s     = lambda_surface_stable
952       lambda_surface_u     = lambda_surface_unstable
953       f_sw_in              = f_shortwave_incoming
954       z0                   = z0_eb
955       z0h                  = z0h_eb
[1788]956       z0q                  = z0q_eb
[1551]957
958!
[1496]959!--    Possibly do user-defined actions (e.g. define heterogeneous land surface)
960       CALL user_init_land_surface
961
[1788]962!
963!--    Set flag parameter if vegetation type was set to a water surface. Also
964!--    set temperature to a constant value in all "soil" layers.
965       DO  i = nxlg, nxrg
966          DO  j = nysg, nyng
967             IF ( veg_type_2d(j,i) == 14  .OR.  veg_type_2d(j,i) == 15 )  THEN
968                water_surface(j,i) = .TRUE.
969             ELSEIF ( veg_type_2d(j,i) == 20 )  THEN
970                pave_surface(j,i) = .TRUE.
971                m_soil(:,j,i) = 0.0_wp
972             ENDIF
[1500]973
[1788]974          ENDDO
975       ENDDO
976
977!
978!--    Calculate new roughness lengths (for water surfaces only)
979       CALL calc_z0_water_surface
980
[1551]981       t_soil_p    = t_soil
982       m_soil_p    = m_soil
983       m_liq_eb_p  = m_liq_eb
[1709]984       t_surface_p = t_surface
[1496]985
[1709]986
987
[1551]988!--    Store initial profiles of t_soil and m_soil (assuming they are
989!--    horizontally homogeneous on this PE)
990       hom(nzb_soil:nzt_soil,1,90,:)  = SPREAD( t_soil(nzb_soil:nzt_soil,      &
991                                                nysg,nxlg), 2,                 &
992                                                statistic_regions+1 )
993       hom(nzb_soil:nzt_soil,1,92,:)  = SPREAD( m_soil(nzb_soil:nzt_soil,      &
994                                                nysg,nxlg), 2,                 &
995                                                statistic_regions+1 )
996
[1496]997    END SUBROUTINE init_lsm
998
999
1000
1001!------------------------------------------------------------------------------!
1002! Description:
1003! ------------
[1682]1004!> Solver for the energy balance at the surface.
[1496]1005!------------------------------------------------------------------------------!
1006    SUBROUTINE lsm_energy_balance
1007
1008
1009       IMPLICIT NONE
1010
[1682]1011       INTEGER(iwp) ::  i         !< running index
1012       INTEGER(iwp) ::  j         !< running index
1013       INTEGER(iwp) ::  k, ks     !< running index
[1496]1014
[1788]1015       REAL(wp) :: c_surface_tmp,& !< temporary variable for storing the volumetric heat capacity of the surface
1016                   f1,          & !< resistance correction term 1
[1682]1017                   f2,          & !< resistance correction term 2
1018                   f3,          & !< resistance correction term 3
1019                   m_min,       & !< minimum soil moisture
1020                   e,           & !< water vapour pressure
1021                   e_s,         & !< water vapour saturation pressure
1022                   e_s_dt,      & !< derivate of e_s with respect to T
1023                   tend,        & !< tendency
1024                   dq_s_dt,     & !< derivate of q_s with respect to T
1025                   coef_1,      & !< coef. for prognostic equation
1026                   coef_2,      & !< coef. for prognostic equation
1027                   f_qsws,      & !< factor for qsws_eb
1028                   f_qsws_veg,  & !< factor for qsws_veg_eb
1029                   f_qsws_soil, & !< factor for qsws_soil_eb
1030                   f_qsws_liq,  & !< factor for qsws_liq_eb
1031                   f_shf,       & !< factor for shf_eb
1032                   lambda_surface, & !< Current value of lambda_surface
[1691]1033                   m_liq_eb_max,   & !< maxmimum value of the liq. water reservoir
[1709]1034                   pt1,         & !< potential temperature at first grid level
1035                   qv1            !< specific humidity at first grid level
[1496]1036
1037!
1038!--    Calculate the exner function for the current time step
1039       exn = ( surface_pressure / 1000.0_wp )**0.286_wp
1040
[1691]1041       DO  i = nxlg, nxrg
1042          DO  j = nysg, nyng
[1500]1043             k = nzb_s_inner(j,i)
[1496]1044
1045!
[1788]1046!--          Set lambda_surface according to stratification between skin layer and soil
1047             IF (  .NOT.  pave_surface(j,i) )  THEN
1048
1049                c_surface_tmp = c_surface
1050
1051                IF ( t_surface(j,i) >= t_soil(nzb_soil,j,i))  THEN
1052                   lambda_surface = lambda_surface_s(j,i)
1053                ELSE
1054                   lambda_surface = lambda_surface_u(j,i)
1055                ENDIF
[1496]1056             ELSE
[1788]1057
1058                c_surface_tmp = pave_heat_capacity * dz_soil(nzb_soil) * 0.5_wp
1059                lambda_surface = pave_heat_conductivity * ddz_soil(nzb_soil)
1060
[1496]1061             ENDIF
[1500]1062
[1496]1063!
[1500]1064!--          First step: calculate aerodyamic resistance. As pt, us, ts
1065!--          are not available for the prognostic time step, data from the last
1066!--          time step is used here. Note that this formulation is the
1067!--          equivalent to the ECMWF formulation using drag coefficients
[1691]1068             IF ( cloud_physics )  THEN
[1709]1069                pt1 = pt(k+1,j,i) + l_d_cp * pt_d_t(k+1) * ql(k+1,j,i)
1070                qv1 = q(k+1,j,i) - ql(k+1,j,i)
[1691]1071             ELSE
[1709]1072                pt1 = pt(k+1,j,i)
1073                qv1 = q(k+1,j,i)
[1691]1074             ENDIF
[1496]1075
[1709]1076             r_a(j,i) = (pt1 - pt(k,j,i)) / (ts(j,i) * us(j,i) + 1.0E-20_wp)
[1691]1077
[1496]1078!
[1709]1079!--          Make sure that the resistance does not drop to zero
1080             IF ( ABS(r_a(j,i)) < 1.0E-10_wp )  r_a(j,i) = 1.0E-10_wp
1081
1082!
[1496]1083!--          Second step: calculate canopy resistance r_canopy
1084!--          f1-f3 here are defined as 1/f1-f3 as in ECMWF documentation
1085 
[1551]1086!--          f1: correction for incoming shortwave radiation (stomata close at
1087!--          night)
[1585]1088             IF ( radiation_scheme /= 'constant' )  THEN
[1709]1089                f1 = MIN( 1.0_wp, ( 0.004_wp * rad_sw_in(k,j,i) + 0.05_wp ) /  &
1090                              (0.81_wp * (0.004_wp * rad_sw_in(k,j,i)          &
[1585]1091                               + 1.0_wp)) )
[1551]1092             ELSE
1093                f1 = 1.0_wp
1094             ENDIF
[1496]1095
[1709]1096
[1496]1097!
[1551]1098!--          f2: correction for soil moisture availability to plants (the
1099!--          integrated soil moisture must thus be considered here)
1100!--          f2 = 0 for very dry soils
[1496]1101             m_total = 0.0_wp
[1691]1102             DO  ks = nzb_soil, nzt_soil
[1551]1103                 m_total = m_total + root_fr(ks,j,i)                           &
1104                           * MAX(m_soil(ks,j,i),m_wilt(j,i))
[1496]1105             ENDDO 
1106
[1788]1107             IF ( m_total > m_wilt(j,i)  .AND.  m_total < m_fc(j,i) )  THEN
[1496]1108                f2 = ( m_total - m_wilt(j,i) ) / (m_fc(j,i) - m_wilt(j,i) )
[1691]1109             ELSEIF ( m_total >= m_fc(j,i) )  THEN
[1551]1110                f2 = 1.0_wp
[1496]1111             ELSE
1112                f2 = 1.0E-20_wp
1113             ENDIF
1114
1115!
1116!--          Calculate water vapour pressure at saturation
[1551]1117             e_s = 0.01_wp * 610.78_wp * EXP( 17.269_wp * ( t_surface(j,i)     &
1118                           - 273.16_wp ) / ( t_surface(j,i) - 35.86_wp ) )
[1496]1119
1120!
1121!--          f3: correction for vapour pressure deficit
[1551]1122             IF ( g_d(j,i) /= 0.0_wp )  THEN
[1496]1123!
1124!--             Calculate vapour pressure
[1709]1125                e  = qv1 * surface_pressure / 0.622_wp
[1551]1126                f3 = EXP ( -g_d(j,i) * (e_s - e) )
[1496]1127             ELSE
1128                f3 = 1.0_wp
1129             ENDIF
1130
1131!
[1551]1132!--          Calculate canopy resistance. In case that c_veg is 0 (bare soils),
1133!--          this calculation is obsolete, as r_canopy is not used below.
[1496]1134!--          To do: check for very dry soil -> r_canopy goes to infinity
[1551]1135             r_canopy(j,i) = r_canopy_min(j,i) / (lai(j,i) * f1 * f2 * f3      &
1136                                             + 1.0E-20_wp)
[1496]1137
1138!
[1551]1139!--          Third step: calculate bare soil resistance r_soil. The Clapp &
1140!--          Hornberger parametrization does not consider c_veg.
[1788]1141             IF ( soil_type_2d(j,i) /= 7 )  THEN
[1551]1142                m_min = c_veg(j,i) * m_wilt(j,i) + (1.0_wp - c_veg(j,i)) *     &
1143                        m_res(j,i)
1144             ELSE
1145                m_min = m_wilt(j,i)
1146             ENDIF
[1496]1147
[1551]1148             f2 = ( m_soil(nzb_soil,j,i) - m_min ) / ( m_fc(j,i) - m_min )
[1513]1149             f2 = MAX(f2,1.0E-20_wp)
[1551]1150             f2 = MIN(f2,1.0_wp)
[1496]1151
1152             r_soil(j,i) = r_soil_min(j,i) / f2
1153
1154!
[1788]1155!--          Calculate the maximum possible liquid water amount on plants and
1156!--          bare surface. For vegetated surfaces, a maximum depth of 0.2 mm is
1157!--          assumed, while paved surfaces might hold up 1 mm of water. The
1158!--          liquid water fraction for paved surfaces is calculated after
1159!--          Noilhan & Planton (1989), while the ECMWF formulation is used for
1160!--          vegetated surfaces and bare soils.
1161             IF ( pave_surface(j,i) )  THEN
1162                m_liq_eb_max = m_max_depth * 5.0_wp
1163                c_liq(j,i) = MIN( 1.0_wp, (m_liq_eb(j,i) / m_liq_eb_max)**0.67 )
1164             ELSE
1165                m_liq_eb_max = m_max_depth * ( c_veg(j,i) * lai(j,i)           &
1166                               + (1.0_wp - c_veg(j,i)) )
1167                c_liq(j,i) = MIN( 1.0_wp, m_liq_eb(j,i) / m_liq_eb_max )
1168             ENDIF
[1496]1169
[1551]1170!
1171!--          Calculate saturation specific humidity
[1496]1172             q_s = 0.622_wp * e_s / surface_pressure
[1500]1173
1174!
[1551]1175!--          In case of dewfall, set evapotranspiration to zero
1176!--          All super-saturated water is then removed from the air
[1788]1177             IF ( humidity  .AND.  q_s <= qv1 )  THEN
[1551]1178                r_canopy(j,i) = 0.0_wp
1179                r_soil(j,i)   = 0.0_wp
[1691]1180             ENDIF
[1496]1181
1182!
1183!--          Calculate coefficients for the total evapotranspiration
[1788]1184!--          In case of water surface, set vegetation and soil fluxes to zero.
1185!--          For pavements, only evaporation of liquid water is possible.
1186             IF ( water_surface(j,i) )  THEN
1187                f_qsws_veg  = 0.0_wp
1188                f_qsws_soil = 0.0_wp
1189                f_qsws_liq  = rho_lv / r_a(j,i)
1190             ELSEIF ( pave_surface (j,i) )  THEN
1191                f_qsws_veg  = 0.0_wp
1192                f_qsws_soil = 0.0_wp
1193                f_qsws_liq  = rho_lv * c_liq(j,i) / r_a(j,i)
1194             ELSE
1195                f_qsws_veg  = rho_lv * c_veg(j,i) * (1.0_wp - c_liq(j,i))/     &
1196                              (r_a(j,i) + r_canopy(j,i))
1197                f_qsws_soil = rho_lv * (1.0_wp - c_veg(j,i)) / (r_a(j,i) +     &
[1551]1198                                                          r_soil(j,i))
[1788]1199                f_qsws_liq  = rho_lv * c_veg(j,i) * c_liq(j,i) / r_a(j,i)
1200             ENDIF
[1500]1201!
1202!--          If soil moisture is below wilting point, plants do no longer
1203!--          transpirate.
[1691]1204!              IF ( m_soil(k,j,i) < m_wilt(j,i) )  THEN
[1551]1205!                 f_qsws_veg = 0.0_wp
1206!              ENDIF
[1496]1207
[1551]1208             f_shf  = rho_cp / r_a(j,i)
1209             f_qsws = f_qsws_veg + f_qsws_soil + f_qsws_liq
1210
1211!
[1496]1212!--          Calculate derivative of q_s for Taylor series expansion
[1571]1213             e_s_dt = e_s * ( 17.269_wp / (t_surface(j,i) - 35.86_wp) -        &
[1551]1214                              17.269_wp*(t_surface(j,i) - 273.16_wp)           &
1215                              / (t_surface(j,i) - 35.86_wp)**2 )
[1496]1216
[1571]1217             dq_s_dt = 0.622_wp * e_s_dt / surface_pressure
[1496]1218
1219!
1220!--          Add LW up so that it can be removed in prognostic equation
[1691]1221             rad_net_l(j,i) = rad_net(j,i) + rad_lw_out(nzb,j,i)
[1496]1222
[1788]1223!
1224!--          Calculate new skin temperature
[1500]1225             IF ( humidity )  THEN
[1691]1226#if defined ( __rrtmg )
[1496]1227!
[1500]1228!--             Numerator of the prognostic equation
[1709]1229                coef_1 = rad_net_l(j,i) + rad_lw_out_change_0(j,i)             &
[1691]1230                         * t_surface(j,i) - rad_lw_out(nzb,j,i)                &
[1709]1231                         + f_shf * pt1 + f_qsws * ( qv1 - q_s                  &
[1571]1232                         + dq_s_dt * t_surface(j,i) ) + lambda_surface         &
[1551]1233                         * t_soil(nzb_soil,j,i)
[1496]1234
1235!
[1500]1236!--             Denominator of the prognostic equation
[1709]1237                coef_2 = rad_lw_out_change_0(j,i) + f_qsws * dq_s_dt           &
[1691]1238                         + lambda_surface + f_shf / exn
1239#else
[1496]1240
[1691]1241!
1242!--             Numerator of the prognostic equation
1243                coef_1 = rad_net_l(j,i) + 3.0_wp * sigma_sb                    &
1244                         * t_surface(j,i) ** 4                                 &
[1709]1245                         + f_shf * pt1 + f_qsws * ( qv1                        &
[1691]1246                         - q_s + dq_s_dt * t_surface(j,i) )                    &
1247                         + lambda_surface * t_soil(nzb_soil,j,i)
1248
1249!
[1709]1250!--             Denominator of the prognostic equation
1251                coef_2 = 4.0_wp * sigma_sb * t_surface(j,i) ** 3 + f_qsws      &
1252                         * dq_s_dt + lambda_surface + f_shf / exn
[1691]1253 
1254#endif
[1500]1255             ELSE
1256
[1691]1257#if defined ( __rrtmg )
[1500]1258!
1259!--             Numerator of the prognostic equation
[1709]1260                coef_1 = rad_net_l(j,i) + rad_lw_out_change_0(j,i)             &
[1691]1261                         * t_surface(j,i) - rad_lw_out(nzb,j,i)                &
[1709]1262                         + f_shf * pt1  + lambda_surface                       &
[1551]1263                         * t_soil(nzb_soil,j,i)
[1500]1264
1265!
1266!--             Denominator of the prognostic equation
[1709]1267                coef_2 = rad_lw_out_change_0(j,i) + lambda_surface + f_shf / exn
[1691]1268#else
1269
1270!
1271!--             Numerator of the prognostic equation
1272                coef_1 = rad_net_l(j,i) + 3.0_wp * sigma_sb                    &
[1709]1273                          * t_surface(j,i) ** 4 + f_shf * pt1                  &
[1691]1274                          + lambda_surface * t_soil(nzb_soil,j,i)
1275
1276!
1277!--             Denominator of the prognostic equation
[1571]1278                coef_2 = 4.0_wp * sigma_sb * t_surface(j,i) ** 3               &
[1551]1279                         + lambda_surface + f_shf / exn
[1691]1280#endif
[1500]1281             ENDIF
1282
[1496]1283             tend = 0.0_wp
1284
1285!
[1551]1286!--          Implicit solution when the surface layer has no heat capacity,
[1496]1287!--          otherwise use RK3 scheme.
[1788]1288             t_surface_p(j,i) = ( coef_1 * dt_3d * tsc(2) + c_surface_tmp *    &
1289                                t_surface(j,i) ) / ( c_surface_tmp + coef_2    &
1290                                * dt_3d * tsc(2) ) 
1291
[1496]1292!
1293!--          Add RK3 term
[1788]1294             IF ( c_surface_tmp /= 0.0_wp )  THEN
[1496]1295
[1691]1296                t_surface_p(j,i) = t_surface_p(j,i) + dt_3d * tsc(3)           &
1297                                   * tt_surface_m(j,i)
1298
[1496]1299!
[1691]1300!--             Calculate true tendency
1301                tend = (t_surface_p(j,i) - t_surface(j,i) - dt_3d * tsc(3)     &
1302                       * tt_surface_m(j,i)) / (dt_3d  * tsc(2))
[1496]1303!
[1691]1304!--             Calculate t_surface tendencies for the next Runge-Kutta step
1305                IF ( timestep_scheme(1:5) == 'runge' )  THEN
1306                   IF ( intermediate_timestep_count == 1 )  THEN
1307                      tt_surface_m(j,i) = tend
1308                   ELSEIF ( intermediate_timestep_count <                      &
1309                            intermediate_timestep_count_max )  THEN
1310                      tt_surface_m(j,i) = -9.5625_wp * tend + 5.3125_wp        &
1311                                          * tt_surface_m(j,i)
1312                   ENDIF
[1496]1313                ENDIF
1314             ENDIF
1315
1316!
[1757]1317!--          In case of fast changes in the skin temperature, it is possible to
1318!--          update the radiative fluxes independently from the prescribed
1319!--          radiation call frequency. This effectively prevents oscillations,
1320!--          especially when setting skip_time_do_radiation /= 0. The threshold
1321!--          value of 0.2 used here is just a first guess. This method should be
1322!--          revised in the future as tests have shown that the threshold is
1323!--          often reached, when no oscillations would occur (causes immense
1324!--          computing time for the radiation code).
[1788]1325             IF ( ABS( t_surface_p(j,i) - t_surface(j,i) ) > 0.2_wp  .AND.     &
[1757]1326                  unscheduled_radiation_calls )  THEN
[1691]1327                force_radiation_call_l = .TRUE.
1328             ENDIF
1329
1330             pt(k,j,i) = t_surface_p(j,i) / exn
1331
1332!
[1496]1333!--          Calculate fluxes
[1691]1334#if defined ( __rrtmg )
[1709]1335             rad_net_l(j,i)   = rad_net_l(j,i) + rad_lw_out_change_0(j,i)      &
[1691]1336                                * t_surface(j,i) - rad_lw_out(nzb,j,i)         &
[1709]1337                                - rad_lw_out_change_0(j,i) * t_surface_p(j,i)
[1691]1338
1339             IF ( rrtm_idrv == 1 )  THEN
1340                rad_net(j,i) = rad_net_l(j,i)
1341                rad_lw_out(nzb,j,i) = rad_lw_out(nzb,j,i)                      &
[1709]1342                                      + rad_lw_out_change_0(j,i)               &
[1691]1343                                      * ( t_surface_p(j,i) - t_surface(j,i) )
1344             ENDIF
1345#else
1346             rad_net_l(j,i)   = rad_net_l(j,i) + 3.0_wp * sigma_sb             &
1347                                * t_surface(j,i)**4 - 4.0_wp * sigma_sb        &
1348                                * t_surface(j,i)**3 * t_surface_p(j,i)
1349#endif
1350
[1551]1351             ghf_eb(j,i)    = lambda_surface * (t_surface_p(j,i)               &
1352                              - t_soil(nzb_soil,j,i))
[1496]1353
[1709]1354             shf_eb(j,i)    = - f_shf * ( pt1 - pt(k,j,i) )
[1691]1355
1356             shf(j,i) = shf_eb(j,i) / rho_cp
1357
[1500]1358             IF ( humidity )  THEN
[1709]1359                qsws_eb(j,i)  = - f_qsws    * ( qv1 - q_s + dq_s_dt            &
[1571]1360                                * t_surface(j,i) - dq_s_dt * t_surface_p(j,i) )
[1496]1361
[1691]1362                qsws(j,i) = qsws_eb(j,i) / rho_lv
1363
[1709]1364                qsws_veg_eb(j,i)  = - f_qsws_veg  * ( qv1 - q_s                &
[1571]1365                                    + dq_s_dt * t_surface(j,i) - dq_s_dt       &
[1551]1366                                    * t_surface_p(j,i) )
1367
[1709]1368                qsws_soil_eb(j,i) = - f_qsws_soil * ( qv1 - q_s                &
[1571]1369                                    + dq_s_dt * t_surface(j,i) - dq_s_dt       &
[1551]1370                                    * t_surface_p(j,i) )
1371
[1709]1372                qsws_liq_eb(j,i)  = - f_qsws_liq  * ( qv1 - q_s                &
[1571]1373                                    + dq_s_dt * t_surface(j,i) - dq_s_dt       &
[1551]1374                                    * t_surface_p(j,i) )
[1500]1375             ENDIF
[1788]1376
[1500]1377!
1378!--          Calculate the true surface resistance
[1691]1379             IF ( qsws_eb(j,i) == 0.0_wp )  THEN
[1551]1380                r_s(j,i) = 1.0E10_wp
[1496]1381             ELSE
[1788]1382                r_s(j,i) = - rho_lv * ( qv1 - q_s + dq_s_dt                    &
[1571]1383                           * t_surface(j,i) - dq_s_dt * t_surface_p(j,i) )     &
[1551]1384                           / qsws_eb(j,i) - r_a(j,i)
[1496]1385             ENDIF
1386
1387!
1388!--          Calculate change in liquid water reservoir due to dew fall or
[1500]1389!--          evaporation of liquid water
1390             IF ( humidity )  THEN
[1496]1391!
[1571]1392!--             If precipitation is activated, add rain water to qsws_liq_eb
1393!--             and qsws_soil_eb according the the vegetation coverage.
[1500]1394!--             precipitation_rate is given in mm.
1395                IF ( precipitation )  THEN
[1571]1396
1397!
1398!--                Add precipitation to liquid water reservoir, if possible.
[1788]1399!--                Otherwise, add the water to soil. In case of
1400!--                pavements, the exceeding water amount is implicitely removed
1401!--                as runoff as qsws_soil_eb is then not used in the soil model
[1691]1402                   IF ( m_liq_eb(j,i) /= m_liq_eb_max )  THEN
[1571]1403                      qsws_liq_eb(j,i) = qsws_liq_eb(j,i)                      &
[1691]1404                                       + c_veg(j,i) * prr(k,j,i) * hyrho(k)    &
[1571]1405                                       * 0.001_wp * rho_l * l_v
1406                   ELSE
1407                      qsws_soil_eb(j,i) = qsws_soil_eb(j,i)                    &
[1691]1408                                        + c_veg(j,i) * prr(k,j,i) * hyrho(k)   &
[1571]1409                                        * 0.001_wp * rho_l * l_v
1410                   ENDIF
1411
1412!--                Add precipitation to bare soil according to the bare soil
1413!--                coverage.
1414                   qsws_soil_eb(j,i) = qsws_soil_eb(j,i) * (1.0_wp             &
[1691]1415                                       - c_veg(j,i)) * prr(k,j,i) * hyrho(k)   &
[1571]1416                                       * 0.001_wp * rho_l * l_v
[1496]1417                ENDIF
[1691]1418
[1500]1419!
1420!--             If the air is saturated, check the reservoir water level
[1691]1421                IF ( qsws_eb(j,i) < 0.0_wp )  THEN
1422
[1500]1423!
[1551]1424!--                Check if reservoir is full (avoid values > m_liq_eb_max)
1425!--                In that case, qsws_liq_eb goes to qsws_soil_eb. In this
1426!--                case qsws_veg_eb is zero anyway (because c_liq = 1),       
1427!--                so that tend is zero and no further check is needed
[1691]1428                   IF ( m_liq_eb(j,i) == m_liq_eb_max )  THEN
[1551]1429                      qsws_soil_eb(j,i) = qsws_soil_eb(j,i)                    &
1430                                           + qsws_liq_eb(j,i)
1431                      qsws_liq_eb(j,i)  = 0.0_wp
[1500]1432                   ENDIF
[1496]1433
1434!
[1551]1435!--                In case qsws_veg_eb becomes negative (unphysical behavior),
1436!--                let the water enter the liquid water reservoir as dew on the
[1500]1437!--                plant
[1691]1438                   IF ( qsws_veg_eb(j,i) < 0.0_wp )  THEN
[1551]1439                      qsws_liq_eb(j,i) = qsws_liq_eb(j,i) + qsws_veg_eb(j,i)
1440                      qsws_veg_eb(j,i) = 0.0_wp
[1500]1441                   ENDIF
1442                ENDIF                   
[1496]1443 
[1551]1444                tend = - qsws_liq_eb(j,i) * drho_l_lv
[1496]1445
[1551]1446                m_liq_eb_p(j,i) = m_liq_eb(j,i) + dt_3d * ( tsc(2) * tend      &
1447                                                   + tsc(3) * tm_liq_eb_m(j,i) )
[1496]1448
1449!
[1500]1450!--             Check if reservoir is overfull -> reduce to maximum
1451!--             (conservation of water is violated here)
[1551]1452                m_liq_eb_p(j,i) = MIN(m_liq_eb_p(j,i),m_liq_eb_max)
[1496]1453
1454!
[1500]1455!--             Check if reservoir is empty (avoid values < 0.0)
1456!--             (conservation of water is violated here)
[1551]1457                m_liq_eb_p(j,i) = MAX(m_liq_eb_p(j,i),0.0_wp)
[1496]1458
1459
1460!
[1551]1461!--             Calculate m_liq_eb tendencies for the next Runge-Kutta step
[1500]1462                IF ( timestep_scheme(1:5) == 'runge' )  THEN
1463                   IF ( intermediate_timestep_count == 1 )  THEN
[1551]1464                      tm_liq_eb_m(j,i) = tend
[1500]1465                   ELSEIF ( intermediate_timestep_count <                      &
1466                            intermediate_timestep_count_max )  THEN
[1551]1467                      tm_liq_eb_m(j,i) = -9.5625_wp * tend + 5.3125_wp         &
1468                                      * tm_liq_eb_m(j,i)
[1500]1469                   ENDIF
[1496]1470                ENDIF
1471
[1500]1472             ENDIF
1473
[1496]1474          ENDDO
[1500]1475       ENDDO
[1496]1476
[1691]1477!
1478!--    Make a logical OR for all processes. Force radiation call if at
[1757]1479!--    least one processor reached the threshold change in skin temperature
[1788]1480       IF ( unscheduled_radiation_calls  .AND.  intermediate_timestep_count    &
[1757]1481            == intermediate_timestep_count_max-1 )  THEN
[1697]1482#if defined( __parallel )
[1691]1483          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1484          CALL MPI_ALLREDUCE( force_radiation_call_l, force_radiation_call,    &
1485                              1, MPI_LOGICAL, MPI_LOR, comm2d, ierr )
1486#else
1487          force_radiation_call = force_radiation_call_l
1488#endif
1489          force_radiation_call_l = .FALSE.
1490       ENDIF
1491
[1788]1492!
1493!--    Calculate surface specific humidity
1494       IF ( humidity )  THEN
1495          CALL calc_q_surface
1496       ENDIF
[1691]1497
[1788]1498!
1499!--    Calculate new roughness lengths (for water surfaces only)
1500       CALL calc_z0_water_surface
1501
1502
[1496]1503    END SUBROUTINE lsm_energy_balance
1504
1505
1506!------------------------------------------------------------------------------!
1507! Description:
1508! ------------
[1682]1509!> Soil model as part of the land surface model. The model predicts soil
1510!> temperature and water content.
[1496]1511!------------------------------------------------------------------------------!
1512    SUBROUTINE lsm_soil_model
1513
1514
1515       IMPLICIT NONE
1516
[1682]1517       INTEGER(iwp) ::  i   !< running index
1518       INTEGER(iwp) ::  j   !< running index
1519       INTEGER(iwp) ::  k   !< running index
[1496]1520
[1682]1521       REAL(wp)     :: h_vg !< Van Genuchten coef. h
[1496]1522
[1682]1523       REAL(wp), DIMENSION(nzb_soil:nzt_soil) :: gamma_temp,  & !< temp. gamma
[1691]1524                                                 lambda_temp, & !< temp. lambda
1525                                                 tend           !< tendency
[1496]1526
[1691]1527       DO  i = nxlg, nxrg
1528          DO  j = nysg, nyng
[1788]1529
1530             IF ( pave_surface(j,i) )  THEN
1531                rho_c_total(nzb_soil,j,i) = pave_heat_capacity
1532                lambda_temp(nzb_soil)     = pave_heat_conductivity
1533             ENDIF
1534
1535             IF (  .NOT.  water_surface(j,i) )  THEN
1536                DO  k = nzb_soil, nzt_soil
1537
1538
1539                   IF ( pave_surface(j,i)  .AND.  zs(k) <= pave_depth )  THEN
1540                   
1541                      rho_c_total(k,j,i) = pave_heat_capacity
1542                      lambda_temp(k)     = pave_heat_conductivity   
1543
1544                   ELSE           
[1496]1545!
[1788]1546!--                   Calculate volumetric heat capacity of the soil, taking
1547!--                   into account water content
1548                      rho_c_total(k,j,i) = (rho_c_soil * (1.0_wp - m_sat(j,i)) &
1549                                           + rho_c_water * m_soil(k,j,i))
[1496]1550
1551!
[1788]1552!--                   Calculate soil heat conductivity at the center of the soil
1553!--                   layers
1554                      lambda_h_sat = lambda_h_sm ** (1.0_wp - m_sat(j,i)) *    &
1555                                     lambda_h_water ** m_soil(k,j,i)
[1691]1556
[1788]1557                      ke = 1.0_wp + LOG10(MAX(0.1_wp,m_soil(k,j,i) / m_sat(j,i)))
[1691]1558
[1788]1559                      lambda_temp(k) = ke * (lambda_h_sat - lambda_h_dry) +    &
1560                                       lambda_h_dry
1561                   ENDIF
[1496]1562
[1788]1563                ENDDO
[1496]1564
1565!
[1788]1566!--             Calculate soil heat conductivity (lambda_h) at the _stag level
1567!--             using linear interpolation. For pavement surface, the
1568!--             true pavement depth is considered
1569                DO  k = nzb_soil, nzt_soil-1
1570                   IF ( pave_surface(j,i)  .AND.  zs(k)   < pave_depth         &
1571                                           .AND.  zs(k+1) > pave_depth )  THEN
1572                      lambda_h(k,j,i) = ( pave_depth - zs(k) ) / dz_soil(k+1)  &
1573                                        * lambda_temp(k)                       &
1574                                        + ( 1.0_wp - ( pave_depth - zs(k) )    &
1575                                        / dz_soil(k+1) ) * lambda_temp(k+1)
1576                   ELSE
1577                      lambda_h(k,j,i) = ( lambda_temp(k+1) + lambda_temp(k) )  &
1578                                        * 0.5_wp
1579                   ENDIF
1580                ENDDO
1581                lambda_h(nzt_soil,j,i) = lambda_temp(nzt_soil)
[1496]1582
[1788]1583
1584
1585
[1496]1586!
[1788]1587!--             Prognostic equation for soil temperature t_soil
1588                tend(:) = 0.0_wp
[1496]1589
[1788]1590                tend(nzb_soil) = (1.0_wp/rho_c_total(nzb_soil,j,i)) *          &
1591                          ( lambda_h(nzb_soil,j,i) * ( t_soil(nzb_soil+1,j,i)  &
1592                            - t_soil(nzb_soil,j,i) ) * ddz_soil(nzb_soil+1)    &
1593                            + ghf_eb(j,i) ) * ddz_soil_stag(nzb_soil)
[1691]1594
[1788]1595                DO  k = nzb_soil+1, nzt_soil
1596                   tend(k) = (1.0_wp/rho_c_total(k,j,i))                       &
1597                             * (   lambda_h(k,j,i)                             &
1598                                 * ( t_soil(k+1,j,i) - t_soil(k,j,i) )         &
1599                                 * ddz_soil(k+1)                               &
1600                                 - lambda_h(k-1,j,i)                           &
1601                                 * ( t_soil(k,j,i) - t_soil(k-1,j,i) )         &
1602                                 * ddz_soil(k)                                 &
1603                               ) * ddz_soil_stag(k)
[1691]1604
[1788]1605                ENDDO
[1691]1606
[1788]1607                t_soil_p(nzb_soil:nzt_soil,j,i) = t_soil(nzb_soil:nzt_soil,j,i)&
1608                                                  + dt_3d * ( tsc(2)           &
1609                                                  * tend(nzb_soil:nzt_soil)    & 
1610                                                  + tsc(3)                     &
1611                                                  * tt_soil_m(:,j,i) )   
[1496]1612
1613!
[1788]1614!--             Calculate t_soil tendencies for the next Runge-Kutta step
1615                IF ( timestep_scheme(1:5) == 'runge' )  THEN
1616                   IF ( intermediate_timestep_count == 1 )  THEN
1617                      DO  k = nzb_soil, nzt_soil
1618                         tt_soil_m(k,j,i) = tend(k)
1619                      ENDDO
1620                   ELSEIF ( intermediate_timestep_count <                      &
1621                            intermediate_timestep_count_max )  THEN
1622                      DO  k = nzb_soil, nzt_soil
1623                         tt_soil_m(k,j,i) = -9.5625_wp * tend(k) + 5.3125_wp   &
[1551]1624                                         * tt_soil_m(k,j,i)
[1788]1625                      ENDDO
1626                   ENDIF
[1496]1627                ENDIF
1628
1629
[1788]1630                DO  k = nzb_soil, nzt_soil
[1551]1631
[1496]1632!
[1788]1633!--                Calculate soil diffusivity at the center of the soil layers
1634                   lambda_temp(k) = (- b_ch * gamma_w_sat(j,i) * psi_sat       &
1635                                     / m_sat(j,i) ) * ( MAX( m_soil(k,j,i),    &
1636                                     m_wilt(j,i) ) / m_sat(j,i) )**(           &
1637                                     b_ch + 2.0_wp )
[1496]1638
1639!
[1788]1640!--                Parametrization of Van Genuchten
1641                   IF ( soil_type /= 7 )  THEN
[1551]1642!
[1788]1643!--                   Calculate the hydraulic conductivity after Van Genuchten
1644!--                   (1980)
1645                      h_vg = ( ( (m_res(j,i) - m_sat(j,i)) / ( m_res(j,i) -    &
1646                                 MAX( m_soil(k,j,i), m_wilt(j,i) ) ) )**(      &
1647                                 n_vg(j,i) / (n_vg(j,i) - 1.0_wp ) ) - 1.0_wp  &
1648                             )**( 1.0_wp / n_vg(j,i) ) / alpha_vg(j,i)
[1496]1649
[1691]1650
[1788]1651                      gamma_temp(k) = gamma_w_sat(j,i) * ( ( (1.0_wp +         &
1652                                      ( alpha_vg(j,i) * h_vg )**n_vg(j,i))**(  &
1653                                      1.0_wp - 1.0_wp / n_vg(j,i) ) - (        &
1654                                      alpha_vg(j,i) * h_vg )**( n_vg(j,i)      &
1655                                      - 1.0_wp) )**2 )                         &
1656                                      / ( ( 1.0_wp + ( alpha_vg(j,i) * h_vg    &
1657                                      )**n_vg(j,i) )**( ( 1.0_wp  - 1.0_wp     &
1658                                      / n_vg(j,i) ) *( l_vg(j,i) + 2.0_wp) ) )
[1496]1659
[1551]1660!
[1788]1661!--                Parametrization of Clapp & Hornberger
1662                   ELSE
1663                      gamma_temp(k) = gamma_w_sat(j,i) * ( m_soil(k,j,i)       &
1664                                      / m_sat(j,i) )**(2.0_wp * b_ch + 3.0_wp)
1665                   ENDIF
[1551]1666
[1788]1667                ENDDO
[1496]1668
1669!
[1788]1670!--             Prognostic equation for soil moisture content. Only performed,
1671!--             when humidity is enabled in the atmosphere and the surface type
1672!--             is not pavement (implies dry soil below).
1673                IF ( humidity  .AND.  .NOT.  pave_surface(j,i) )  THEN
1674!
1675!--                Calculate soil diffusivity (lambda_w) at the _stag level
1676!--                using linear interpolation. To do: replace this with
1677!--                ECMWF-IFS Eq. 8.81
1678                   DO  k = nzb_soil, nzt_soil-1
[1496]1679                     
[1788]1680                      lambda_w(k,j,i) = ( lambda_temp(k+1) + lambda_temp(k) )  &
1681                                        * 0.5_wp
1682                      gamma_w(k,j,i)  = ( gamma_temp(k+1) + gamma_temp(k) )    &
1683                                        * 0.5_wp
[1496]1684
[1788]1685                   ENDDO
[1496]1686
1687!
1688!
[1788]1689!--                In case of a closed bottom (= water content is conserved), set
1690!--                hydraulic conductivity to zero to that no water will be lost
1691!--                in the bottom layer.
1692                   IF ( conserve_water_content )  THEN
1693                      gamma_w(nzt_soil,j,i) = 0.0_wp
1694                   ELSE
1695                      gamma_w(nzt_soil,j,i) = gamma_temp(nzt_soil)
1696                   ENDIF     
[1496]1697
[1788]1698!--                The root extraction (= root_extr * qsws_veg_eb / (rho_l * l_v))
1699!--                ensures the mass conservation for water. The transpiration of
1700!--                plants equals the cumulative withdrawals by the roots in the
1701!--                soil. The scheme takes into account the availability of water
1702!--                in the soil layers as well as the root fraction in the
1703!--                respective layer. Layer with moisture below wilting point will
1704!--                not contribute, which reflects the preference of plants to
1705!--                take water from moister layers.
[1496]1706
1707!
[1788]1708!--                Calculate the root extraction (ECMWF 7.69, the sum of root_extr
1709!--                = 1). The energy balance solver guarantees a positive
1710!--                transpiration, so that there is no need for an additional check.
[1691]1711                   DO  k = nzb_soil, nzt_soil
[1788]1712                       IF ( m_soil(k,j,i) > m_wilt(j,i) )  THEN
1713                          m_total = m_total + root_fr(k,j,i) * m_soil(k,j,i)
1714                       ENDIF
1715                   ENDDO 
[1496]1716
[1788]1717                   IF ( m_total > 0.0_wp )  THEN
1718                      DO  k = nzb_soil, nzt_soil
1719                         IF ( m_soil(k,j,i) > m_wilt(j,i) )  THEN
1720                            root_extr(k) = root_fr(k,j,i) * m_soil(k,j,i) / m_total
1721                         ELSE
1722                            root_extr(k) = 0.0_wp
1723                         ENDIF
1724                      ENDDO
1725                   ENDIF
1726
[1496]1727!
[1788]1728!--                Prognostic equation for soil water content m_soil.
1729                   tend(:) = 0.0_wp
1730
1731                   tend(nzb_soil) = ( lambda_w(nzb_soil,j,i) * (               &
[1551]1732                            m_soil(nzb_soil+1,j,i) - m_soil(nzb_soil,j,i) )    &
[1691]1733                            * ddz_soil(nzb_soil+1) - gamma_w(nzb_soil,j,i) - ( &
[1788]1734                               root_extr(nzb_soil) * qsws_veg_eb(j,i)          &
1735                               + qsws_soil_eb(j,i) ) * drho_l_lv )             &
1736                               * ddz_soil_stag(nzb_soil)
[1496]1737
[1788]1738                   DO  k = nzb_soil+1, nzt_soil-1
1739                      tend(k) = ( lambda_w(k,j,i) * ( m_soil(k+1,j,i)          &
1740                                - m_soil(k,j,i) ) * ddz_soil(k+1)              &
1741                                - gamma_w(k,j,i)                               &
1742                                - lambda_w(k-1,j,i) * (m_soil(k,j,i) -         &
1743                                m_soil(k-1,j,i)) * ddz_soil(k)                 &
1744                                + gamma_w(k-1,j,i) - (root_extr(k)             &
1745                                * qsws_veg_eb(j,i) * drho_l_lv)                &
1746                                ) * ddz_soil_stag(k)
[1496]1747
[1788]1748                   ENDDO
1749                   tend(nzt_soil) = ( - gamma_w(nzt_soil,j,i)                  &
1750                                           - lambda_w(nzt_soil-1,j,i)          &
1751                                           * (m_soil(nzt_soil,j,i)             &
1752                                           - m_soil(nzt_soil-1,j,i))           &
1753                                           * ddz_soil(nzt_soil)                &
1754                                           + gamma_w(nzt_soil-1,j,i) - (       &
1755                                             root_extr(nzt_soil)               &
1756                                           * qsws_veg_eb(j,i) * drho_l_lv  )   &
1757                                     ) * ddz_soil_stag(nzt_soil)             
[1496]1758
[1788]1759                   m_soil_p(nzb_soil:nzt_soil,j,i) = m_soil(nzb_soil:nzt_soil,j,i)&
1760                                                   + dt_3d * ( tsc(2) * tend(:)   &
1761                                                   + tsc(3) * tm_soil_m(:,j,i) )   
1762   
[1496]1763!
[1788]1764!--                Account for dry soils (find a better solution here!)
1765                   DO  k = nzb_soil, nzt_soil
1766                      IF ( m_soil_p(k,j,i) < 0.0_wp )  m_soil_p(k,j,i) = 0.0_wp
1767                   ENDDO
[1496]1768
1769!
[1788]1770!--                Calculate m_soil tendencies for the next Runge-Kutta step
1771                   IF ( timestep_scheme(1:5) == 'runge' )  THEN
1772                      IF ( intermediate_timestep_count == 1 )  THEN
1773                         DO  k = nzb_soil, nzt_soil
1774                            tm_soil_m(k,j,i) = tend(k)
1775                         ENDDO
1776                      ELSEIF ( intermediate_timestep_count <                   &
1777                               intermediate_timestep_count_max )  THEN
1778                         DO  k = nzb_soil, nzt_soil
1779                            tm_soil_m(k,j,i) = -9.5625_wp * tend(k) + 5.3125_wp&
[1496]1780                                     * tm_soil_m(k,j,i)
[1788]1781                         ENDDO
1782                      ENDIF
[1496]1783                   ENDIF
[1788]1784
[1496]1785                ENDIF
1786
1787             ENDIF
1788
1789          ENDDO
1790       ENDDO
1791
[1788]1792    END SUBROUTINE lsm_soil_model
1793
1794
1795!------------------------------------------------------------------------------!
1796! Description:
1797! ------------
1798!> Calculation of roughness length for open water (lakes, ocean). The
1799!> parameterization follows Charnock (1955). Two different implementations
1800!> are available: as in ECMWF-IFS (Beljaars 1994) or as in FLake (Subin et al.
1801!> 2012)
1802!------------------------------------------------------------------------------!
1803    SUBROUTINE calc_z0_water_surface
1804
1805       USE control_parameters,                                                 &
1806           ONLY: g, kappa, molecular_viscosity
1807
1808       IMPLICIT NONE
1809
1810       INTEGER :: i  !< running index
1811       INTEGER :: j  !< running index
1812
1813       REAL(wp), PARAMETER :: alpha_ch  = 0.018_wp !< Charnock constant (0.01-0.11). Use 0.01 for FLake and 0.018 for ECMWF
1814!       REAL(wp), PARAMETER :: pr_number = 0.71_wp !< molecular Prandtl number in the Charnock parameterization (differs from prandtl_number)
1815!       REAL(wp), PARAMETER :: sc_number = 0.66_wp !< molecular Schmidt number in the Charnock parameterization
1816!       REAL(wp) :: re_0 !< near-surface roughness Reynolds number
1817
1818
1819       DO  i = nxlg, nxrg   
1820          DO  j = nysg, nyng
1821             IF ( water_surface(j,i) )  THEN
1822
[1496]1823!
[1788]1824!--             Disabled: FLake parameterization. Ideally, the Charnock
1825!--             coefficient should depend on the water depth and the fetch
1826!--             length
1827!                re_0 = z0(j,i) * us(j,i) / molecular_viscosity
1828!       
1829!                z0(j,i) = MAX( 0.1_wp * molecular_viscosity / us(j,i),            &
1830!                              alpha_ch * us(j,i) / g )
1831!
1832!                z0h(j,i) = z0(j,i) * EXP( - kappa / pr_number * ( 4.0_wp * SQRT( re_0 ) - 3.2_wp ) )
1833!                z0q(j,i) = z0(j,i) * EXP( - kappa / pr_number * ( 4.0_wp * SQRT( re_0 ) - 4.2_wp ) )
[1496]1834
[1788]1835!
1836!--              Set minimum roughness length for u* > 0.2
1837!                IF ( us(j,i) > 0.2_wp )  THEN
1838!                   z0h(j,i) = MAX( 1.0E-5_wp, z0h(j,i) )
1839!                   z0q(j,i) = MAX( 1.0E-5_wp, z0q(j,i) )
1840!                ENDIF
[1496]1841
[1788]1842!
1843!--             ECMWF IFS model parameterization after Beljaars (1994). At low
1844!--             wind speed, the sea surface becomes aerodynamically smooth and
1845!--             the roughness scales with the viscosity. At high wind speed, the
1846!--             Charnock relation is used.
1847                z0(j,i) =   ( 0.11_wp * molecular_viscosity / us(j,i) )        &
1848                          + ( alpha_ch * us(j,i)**2 / g )
[1496]1849
[1788]1850                z0h(j,i) = 0.40_wp * molecular_viscosity / us(j,i)
1851                z0q(j,i) = 0.62_wp * molecular_viscosity / us(j,i)
[1496]1852
[1788]1853             ENDIF
1854          ENDDO
1855       ENDDO
1856
1857    END SUBROUTINE calc_z0_water_surface
1858
1859
[1496]1860!------------------------------------------------------------------------------!
1861! Description:
1862! ------------
[1788]1863!> Calculation of specific humidity of the skin layer (surface). It is assumend
1864!> that the skin is always saturated.
[1496]1865!------------------------------------------------------------------------------!
[1551]1866    SUBROUTINE calc_q_surface
[1496]1867
1868       IMPLICIT NONE
1869
[1682]1870       INTEGER :: i              !< running index
1871       INTEGER :: j              !< running index
1872       INTEGER :: k              !< running index
[1788]1873
[1682]1874       REAL(wp) :: resistance    !< aerodynamic and soil resistance term
[1496]1875
[1691]1876       DO  i = nxlg, nxrg   
1877          DO  j = nysg, nyng
[1496]1878             k = nzb_s_inner(j,i)
1879
1880!
1881!--          Calculate water vapour pressure at saturation
[1691]1882             e_s = 0.01_wp * 610.78_wp * EXP( 17.269_wp * ( t_surface_p(j,i)   &
1883                   - 273.16_wp ) / ( t_surface_p(j,i) - 35.86_wp ) )
[1496]1884
1885!
1886!--          Calculate specific humidity at saturation
1887             q_s = 0.622_wp * e_s / surface_pressure
1888
1889             resistance = r_a(j,i) / (r_a(j,i) + r_s(j,i))
1890
1891!
1892!--          Calculate specific humidity at surface
[1691]1893             IF ( cloud_physics )  THEN
1894                q(k,j,i) = resistance * q_s + (1.0_wp - resistance)            &
1895                             * ( q(k+1,j,i) - ql(k+1,j,i) )
1896             ELSE
1897                q(k,j,i) = resistance * q_s + (1.0_wp - resistance)            &
1898                             * q(k+1,j,i)
1899             ENDIF
[1496]1900
[1691]1901!
1902!--          Update virtual potential temperature
1903             vpt(k,j,i) = pt(k,j,i) * ( 1.0_wp + 0.61_wp * q(k,j,i) )
1904
[1496]1905          ENDDO
1906       ENDDO
1907
[1551]1908    END SUBROUTINE calc_q_surface
[1496]1909
[1788]1910
[1551]1911!------------------------------------------------------------------------------!
1912! Description:
1913! ------------
[1682]1914!> Swapping of timelevels
[1551]1915!------------------------------------------------------------------------------!
1916    SUBROUTINE lsm_swap_timelevel ( mod_count )
[1496]1917
[1551]1918       IMPLICIT NONE
1919
1920       INTEGER, INTENT(IN) :: mod_count
1921
1922#if defined( __nopointer )
1923
1924       t_surface    = t_surface_p
1925       t_soil       = t_soil_p
1926       IF ( humidity )  THEN
1927          m_soil    = m_soil_p
1928          m_liq_eb  = m_liq_eb_p
1929       ENDIF
1930
1931#else
1932   
1933       SELECT CASE ( mod_count )
1934
1935          CASE ( 0 )
1936
[1585]1937             t_surface  => t_surface_1; t_surface_p  => t_surface_2
1938             t_soil     => t_soil_1;    t_soil_p     => t_soil_2
[1551]1939             IF ( humidity )  THEN
[1585]1940                m_soil    => m_soil_1;   m_soil_p    => m_soil_2
1941                m_liq_eb  => m_liq_eb_1; m_liq_eb_p  => m_liq_eb_2
[1551]1942             ENDIF
1943
1944
1945          CASE ( 1 )
1946
[1585]1947             t_surface  => t_surface_2; t_surface_p  => t_surface_1
1948             t_soil     => t_soil_2;    t_soil_p     => t_soil_1
[1551]1949             IF ( humidity )  THEN
[1585]1950                m_soil    => m_soil_2;   m_soil_p    => m_soil_1
1951                m_liq_eb  => m_liq_eb_2; m_liq_eb_p  => m_liq_eb_1
[1551]1952             ENDIF
1953
1954       END SELECT
1955#endif
1956
1957    END SUBROUTINE lsm_swap_timelevel
1958
1959
[1496]1960 END MODULE land_surface_model_mod
Note: See TracBrowser for help on using the repository browser.