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

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

bugfixes in output of land surface variables

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