source: palm/trunk/SOURCE/land_surface_model_mod.f90 @ 1818

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

last commit documented / copyright update

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