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

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

bugfix in land surface / radiation model

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