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

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

further modularization of land surface model

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