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

Last change on this file since 1854 was 1852, checked in by hoffmann, 9 years ago

last commit documented

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