source: palm/trunk/SOURCE/land_surface_model.f90 @ 1553

Last change on this file since 1553 was 1553, checked in by maronga, 10 years ago

further adjustments in land surface model

  • Property svn:keywords set to Id
File size: 68.7 KB
Line 
1 MODULE land_surface_model_mod
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2014 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22! Improved better treatment of roughness lengths. Added default soil temperature
23! profile
24!
25! Former revisions:
26! -----------------
27! $Id: land_surface_model.f90 1553 2015-03-03 17:33:54Z maronga $
28! Bugfix: REAL constants provided with KIND-attribute in call of
29! 1551 2015-03-03 14:18:16Z maronga
30! Flux calculation is now done in prandtl_fluxes. Added support for data output.
31! Vertical indices have been replaced. Restart runs are now possible. Some
32! variables have beem renamed. Bugfix in the prognostic equation for the surface
33! temperature. Introduced z0_eb and z0h_eb, which overwrite the setting of
34! roughness_length and z0_factor. Added Clapp & Hornberger parametrization for
35! the hydraulic conductivity. Bugfix for root fraction and extraction
36! calculation
37!
38! intrinsic function MAX and MIN
39!
40! 1500 2014-12-03 17:42:41Z maronga
41! Corrected calculation of aerodynamic resistance (r_a).
42! Precipitation is now added to liquid water reservoir using LE_liq.
43! Added support for dry runs.
44!
45! 1496 2014-12-02 17:25:50Z maronga
46! Initial revision
47!
48!
49! Description:
50! ------------
51! Land surface model, consisting of a solver for the energy balance at the
52! surface and a four layer soil scheme. The scheme is similar to the TESSEL
53! scheme implemented in the ECMWF IFS model, with modifications according to
54! H-TESSEL. The implementation is based on the formulation implemented in the
55! DALES and UCLA-LES models.
56!
57! To do list:
58! -----------
59! - Check dewfall parametrization for fog simulations
60! - Consider partial absorption of the net shortwave radiation by the surface layer
61! - Allow for water surfaces, check performance for bare soils
62! - Invert indices (running from -3 to 0. Currently: nzb_soil=0, nzt_soil=3))
63! - Implement surface runoff model (required when performing long-term LES with
64!   considerable precipitation
65! Notes:
66! ------
67! - No time step criterion is required as long as the soil layers do not become
68!   too thin
69!------------------------------------------------------------------------------!
70     USE arrays_3d,                                                            &
71         ONLY:  pt, pt_p, q, q_p, qsws, rif, shf, ts, us, z0, z0h
72
73     USE cloud_parameters,                                                     &
74         ONLY:  cp, l_d_r, l_v, precipitation_rate, rho_l, r_d, r_v
75
76     USE control_parameters,                                                   &
77         ONLY:  cloud_physics, dt_3d, humidity, intermediate_timestep_count,   &
78                initializing_actions, intermediate_timestep_count_max,         &
79                max_masks, precipitation, pt_surface, rho_surface,             &
80                roughness_length, surface_pressure, timestep_scheme, tsc,      &
81                z0h_factor
82
83     USE indices,                                                              &
84         ONLY:  nxlg, nxrg, nyng, nysg, nzb_s_inner 
85
86     USE kinds
87
88    USE netcdf_control,                                                        &
89        ONLY:  dots_label, dots_num, dots_unit
90
91     USE radiation_model_mod,                                                  &
92         ONLY:  irad_scheme, rad_net, rad_sw_in, sigma_SB
93
94     USE statistics,                                                           &
95         ONLY:  hom, statistic_regions
96
97    IMPLICIT NONE
98
99!
100!-- LSM model constants
101    INTEGER(iwp), PARAMETER :: nzb_soil = 0, & !: bottom of the soil model (to be switched)
102                               nzt_soil = 3, & !: top of the soil model (to be switched)
103                               nzs = 4         !: number of soil layers (fixed for now)
104
105    INTEGER(iwp) :: dots_soil = 0  !: starting index for timeseries output
106
107    INTEGER(iwp), DIMENSION(0:1) :: id_dim_zs_xy, id_dim_zs_xz, id_dim_zs_yz,  &
108                                    id_dim_zs_3d, id_var_zs_xy,                & 
109                                    id_var_zs_xz, id_var_zs_yz, id_var_zs_3d
110                                   
111    INTEGER(iwp), DIMENSION(1:max_masks,0:1) :: id_dim_zs_mask, id_var_zs_mask
112
113    REAL(wp), PARAMETER ::                     &
114              b_ch               = 6.04_wp,    & ! Clapp & Hornberger exponent
115              lambda_h_dry       = 0.19_wp,    & ! heat conductivity for dry soil   
116              lambda_h_sm        = 3.44_wp,    & ! heat conductivity of the soil matrix
117              lambda_h_water     = 0.57_wp,    & ! heat conductivity of water
118              psi_sat            = -0.388_wp,  & ! soil matrix potential at saturation
119              rho_c_soil         = 2.19E6_wp,  & ! volumetric heat capacity of soil
120              rho_c_water        = 4.20E6_wp,  & ! volumetric heat capacity of water
121              m_max_depth        = 0.0002_wp     ! Maximum capacity of the water reservoir (m)
122
123
124!
125!-- LSM variables
126    INTEGER(iwp) :: veg_type  = 2, & !: vegetation type, 0: user-defined, 1-19: generic (see list)
127                    soil_type = 3    !: soil type, 0: user-defined, 1-6: generic (see list)
128
129    LOGICAL :: conserve_water_content = .TRUE., & !: open or closed bottom surface for the soil model
130               dewfall = .TRUE.,                & !: allow/inhibit dewfall
131               land_surface = .FALSE.             !: flag parameter indicating wheather the lsm is used
132
133!   value 9999999.9_wp -> generic available or user-defined value must be set
134!   otherwise -> no generic variable and user setting is optional
135    REAL(wp) :: alpha_vangenuchten = 9999999.9_wp,      & !: NAMELIST alpha_vg
136                canopy_resistance_coefficient = 9999999.9_wp, & !: NAMELIST g_d
137                c_surface   = 20000.0_wp,               & !: Surface (skin) heat capacity
138                drho_l_lv,                              & !: (rho_l * l_v)**-1
139                exn,                                    & !: value of the Exner function
140                e_s = 0.0_wp,                           & !: saturation water vapour pressure
141                field_capacity = 9999999.9_wp,          & !: NAMELIST m_fc
142                f_shortwave_incoming = 9999999.9_wp,    & !: NAMELIST f_sw_in
143                hydraulic_conductivity = 9999999.9_wp,  & !: NAMELIST gamma_w_sat
144                ke = 0.0_wp,                            & !: Kersten number
145                lambda_surface_stable = 9999999.9_wp,   & !: NAMELIST lambda_surface_s
146                lambda_surface_unstable = 9999999.9_wp, & !: NAMELIST lambda_surface_u
147                leaf_area_index = 9999999.9_wp,         & !: NAMELIST lai
148                l_vangenuchten = 9999999.9_wp,          & !: NAMELIST l_vg
149                min_canopy_resistance = 9999999.9_wp,   & !: NAMELIST r_canopy_min
150                min_soil_resistance = 50.0_wp,          & !: NAMELIST r_soil_min
151                m_total = 0.0_wp,                       & !: weighted total water content of the soil (m3/m3)
152                n_vangenuchten = 9999999.9_wp,          & !: NAMELIST n_vg
153                q_s = 0.0_wp,                           & !: saturation specific humidity
154                residual_moisture = 9999999.9_wp,       & !: NAMELIST m_res
155                rho_cp,                                 & !: rho_surface * cp
156                rho_lv,                                 & !: rho * l_v
157                rd_d_rv,                                & !: r_d / r_v
158                saturation_moisture = 9999999.9_wp,     & !: NAMELIST m_sat
159                vegetation_coverage = 9999999.9_wp,     & !: NAMELIST c_veg
160                wilting_point = 9999999.9_wp,           & !: NAMELIST m_wilt
161                z0_eb  = 9999999.9_wp,                  & !: NAMELIST z0 (lsm_par)
162                z0h_eb = 9999999.9_wp                    !: NAMELIST z0h (lsm_par)
163
164    REAL(wp), DIMENSION(nzb_soil:nzt_soil) :: &
165              ddz_soil,                     & !: 1/dz_soil
166              ddz_soil_stag,                & !: 1/dz_soil_stag
167              dz_soil,                      & !: soil grid spacing (center-center)
168              dz_soil_stag,                 & !: soil grid spacing (edge-edge)
169              root_extr = 0.0_wp,           & !: root extraction
170              root_fraction = (/9999999.9_wp, 9999999.9_wp,    &
171                                9999999.9_wp, 9999999.9_wp /), & !: distribution of root surface area to the individual soil layers
172              zs = (/0.07_wp, 0.28_wp, 1.00_wp,  2.89_wp/),    & !: soil layer depths (m)
173              soil_moisture = 0.0_wp          !: soil moisture content (m3/m3)
174
175    REAL(wp), DIMENSION(nzb_soil:nzt_soil+1) ::   &
176              soil_temperature = (/290.0_wp, 287.0_wp, 285.0_wp,  283.0_wp,    &
177                                   283.0_wp /) !: soil temperature (K)
178
179#if defined( __nopointer )
180    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET :: t_surface,   & !: surface temperature (K)
181                                                     t_surface_p, & !: progn. surface temperature (K)
182                                                     m_liq_eb,    & !: liquid water reservoir (m)
183                                                     m_liq_eb_av, & !: liquid water reservoir (m)
184                                                     m_liq_eb_p     !: progn. liquid water reservoir (m)
185#else
186    REAL(wp), DIMENSION(:,:), POINTER :: t_surface,      &
187                                         t_surface_p,    & 
188                                         m_liq_eb,       & 
189                                         m_liq_eb_p
190
191    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET :: t_surface_1, t_surface_2, &
192                                                     m_liq_eb_av,              &
193                                                     m_liq_eb_1, m_liq_eb_2
194#endif
195
196!
197!-- Temporal tendencies for time stepping           
198    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: tt_surface_m,  & !: surface temperature tendency (K)
199                                             tm_liq_eb_m      !: liquid water reservoir tendency (m)
200
201!
202!-- Energy balance variables               
203    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::                                   &
204              alpha_vg,         & !: coef. of Van Genuchten
205              c_liq,            & !: liquid water coverage (of vegetated area)
206              c_liq_av,         & !: average of c_liq
207              c_soil_av,        & !: average of c_soil
208              c_veg,            & !: vegetation coverage
209              c_veg_av,         & !: average of c_veg
210              f_sw_in,          & !: fraction of absorbed shortwave radiation by the surface layer (not implemented yet)
211              ghf_eb,           & !: ground heat flux
212              ghf_eb_av,        & !: average of ghf_eb
213              gamma_w_sat,      & !: hydraulic conductivity at saturation
214              g_d,              & !: coefficient for dependence of r_canopy on water vapour pressure deficit
215              lai,              & !: leaf area index
216              lai_av,           & !: average of lai
217              lambda_h_sat,     & !: heat conductivity for dry soil
218              lambda_surface_s,    & !: coupling between surface and soil (depends on vegetation type)
219              lambda_surface_u,    & !: coupling between surface and soil (depends on vegetation type)
220              l_vg,             & !: coef. of Van Genuchten
221              m_fc,             & !: soil moisture at field capacity (m3/m3)
222              m_res,            & !: residual soil moisture
223              m_sat,            & !: saturation soil moisture (m3/m3)
224              m_wilt,           & !: soil moisture at permanent wilting point (m3/m3)
225              n_vg,             & !: coef. Van Genuchten 
226              qsws_eb,          & !: surface flux of latent heat (total)
227              qsws_eb_av,       & !: average of qsws_eb
228              qsws_liq_eb,      & !: surface flux of latent heat (liquid water portion)
229              qsws_liq_eb_av,   & !: average of qsws_liq_eb
230              qsws_soil_eb,     & !: surface flux of latent heat (soil portion)
231              qsws_soil_eb_av,  & !: average of qsws_soil_eb
232              qsws_veg_eb,      & !: surface flux of latent heat (vegetation portion)
233              qsws_veg_eb_av,   & !: average of qsws_veg_eb
234              r_a,              & !: aerodynamic resistance
235              r_canopy,         & !: canopy resistance
236              r_soil,           & !: soil resitance
237              r_soil_min,       & !: minimum soil resistance
238              r_s,              & !: total surface resistance (combination of r_soil and r_canopy)         
239              r_canopy_min,     & !: minimum canopy (stomatal) resistance
240              shf_eb,           & !: surface flux of sensible heat
241              shf_eb_av           !: average of shf_eb
242
243    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::                                 &
244              lambda_h, &   !: heat conductivity of soil (?)                           
245              lambda_w, &   !: hydraulic diffusivity of soil (?)
246              gamma_w,  &   !: hydraulic conductivity of soil (?)
247              rho_c_total   !: volumetric heat capacity of the actual soil matrix (?)
248
249#if defined( __nopointer )
250    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                         &
251              t_soil,    & !: Soil temperature (K)
252              t_soil_av, & !: Average of t_soil
253              t_soil_p,  & !: Prog. soil temperature (K)
254              m_soil,    & !: Soil moisture (m3/m3)
255              m_soil_av, & !: Average of m_soil
256              m_soil_p     !: Prog. soil moisture (m3/m3)
257#else
258    REAL(wp), DIMENSION(:,:,:), POINTER ::                                     &
259              t_soil, t_soil_p, &
260              m_soil, m_soil_p   
261
262    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                         &
263              t_soil_av, t_soil_1, t_soil_2,                                   &
264              m_soil_av, m_soil_1, m_soil_2
265
266
267#endif
268
269
270    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::                                 &
271              tt_soil_m, & !: t_soil storage array
272              tm_soil_m, & !: m_soil storage array
273              root_fr      !: root fraction (sum=1)
274
275
276!
277!-- Predefined Land surface classes (veg_type)
278    CHARACTER(26), DIMENSION(0:19) :: veg_type_name = (/          &
279                                   'user defined',                & ! 0
280                                   'crops, mixed farming',        & !  1
281                                   'short grass',                 & !  2
282                                   'evergreen needleleaf trees',  & !  3
283                                   'deciduous needleleaf trees',  & !  4
284                                   'evergreen broadleaf trees' ,  & !  5
285                                   'deciduous broadleaf trees',   & !  6
286                                   'tall grass',                  & !  7
287                                   'desert',                      & !  8
288                                   'tundra',                      & !  9
289                                   'irrigated crops',             & ! 10
290                                   'semidesert',                  & ! 11
291                                   'ice caps and glaciers' ,      & ! 12
292                                   'bogs and marshes',            & ! 13
293                                   'inland water',                & ! 14
294                                   'ocean',                       & ! 15
295                                   'evergreen shrubs',            & ! 16
296                                   'deciduous shrubs',            & ! 17
297                                   'mixed forest/woodland',       & ! 18
298                                   'interrupted forest'           & ! 19
299                                                       /)
300
301!
302!-- Soil model classes (soil_type)
303    CHARACTER(12), DIMENSION(0:7) :: soil_type_name = (/ &
304                                   'user defined',        & ! 0
305                                   'coarse',              & ! 1
306                                   'medium',              & ! 2
307                                   'medium-fine',         & ! 3
308                                   'fine',                & ! 4
309                                   'very fine' ,          & ! 5
310                                   'organic',             & ! 6
311                                   'loamy (CH)'           & ! 7
312                                                        /)
313!
314!-- Land surface parameters according to the respective classes (veg_type)
315
316!
317!-- Land surface parameters I
318!--                          r_canopy_min,     lai,   c_veg,     g_d
319    REAL(wp), DIMENSION(0:3,1:19) :: veg_pars = RESHAPE( (/           &
320                                 180.0_wp, 3.00_wp, 0.90_wp, 0.00_wp, & !  1
321                                 110.0_wp, 2.00_wp, 0.85_wp, 0.00_wp, & !  2
322                                 500.0_wp, 5.00_wp, 0.90_wp, 0.03_wp, & !  3
323                                 500.0_wp, 5.00_wp, 0.90_wp, 0.03_wp, & !  4
324                                 175.0_wp, 5.00_wp, 0.90_wp, 0.03_wp, & !  5
325                                 240.0_wp, 6.00_wp, 0.99_wp, 0.13_wp, & !  6
326                                 100.0_wp, 2.00_wp, 0.70_wp, 0.00_wp, & !  7
327                                 250.0_wp, 0.05_wp, 0.00_wp, 0.00_wp, & !  8
328                                  80.0_wp, 1.00_wp, 0.50_wp, 0.00_wp, & !  9
329                                 180.0_wp, 3.00_wp, 0.90_wp, 0.00_wp, & ! 10
330                                 150.0_wp, 0.50_wp, 0.10_wp, 0.00_wp, & ! 11
331                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 12
332                                 240.0_wp, 4.00_wp, 0.60_wp, 0.00_wp, & ! 13
333                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 14
334                                   0.0_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 15
335                                 225.0_wp, 3.00_wp, 0.50_wp, 0.00_wp, & ! 16
336                                 225.0_wp, 1.50_wp, 0.50_wp, 0.00_wp, & ! 17
337                                 250.0_wp, 5.00_wp, 0.90_wp, 0.03_wp, & ! 18
338                                 175.0_wp, 2.50_wp, 0.90_wp, 0.03_wp  & ! 19
339                                 /), (/ 4, 19 /) )
340
341!
342!-- Land surface parameters II          z0,         z0h
343    REAL(wp), DIMENSION(0:1,1:19) :: roughness_par = RESHAPE( (/ & 
344                                   0.25_wp,  0.25E-2_wp,         & !  1
345                                   0.20_wp,  0.20E-2_wp,         & !  2
346                                   2.00_wp,     2.00_wp,         & !  3
347                                   2.00_wp,     2.00_wp,         & !  4
348                                   2.00_wp,     2.00_wp,         & !  5
349                                   2.00_wp,     2.00_wp,         & !  6
350                                   0.47_wp,  0.47E-2_wp,         & !  7
351                                  0.013_wp, 0.013E-2_wp,         & !  8
352                                  0.034_wp, 0.034E-2_wp,         & !  9
353                                    0.5_wp,  0.50E-2_wp,         & ! 10
354                                   0.17_wp,  0.17E-2_wp,         & ! 11
355                                 1.3E-3_wp,   1.3E-4_wp,         & ! 12
356                                   0.83_wp,  0.83E-2_wp,         & ! 13
357                                   0.00_wp,  0.00E-2_wp,         & ! 14
358                                   0.00_wp,  0.00E-2_wp,         & ! 15
359                                   0.10_wp,  0.10E-2_wp,         & ! 16
360                                   0.25_wp,  0.25E-2_wp,         & ! 17
361                                   2.00_wp,  2.00E-2_wp,         & ! 18
362                                   1.10_wp,  1.10E-2_wp          & ! 19
363                                 /), (/ 2, 19 /) )
364
365!
366!-- Land surface parameters III lambda_surface_s, lambda_surface_u, f_sw_in
367    REAL(wp), DIMENSION(0:2,1:19) :: surface_pars = RESHAPE( (/           &
368                                      10.0_wp,       10.0_wp, 0.05_wp, & !  1
369                                      10.0_wp,       10.0_wp, 0.05_wp, & !  2
370                                      20.0_wp,       15.0_wp, 0.03_wp, & !  3
371                                      20.0_wp,       15.0_wp, 0.03_wp, & !  4
372                                      20.0_wp,       15.0_wp, 0.03_wp, & !  5
373                                      20.0_wp,       15.0_wp, 0.03_wp, & !  6
374                                      10.0_wp,       10.0_wp, 0.05_wp, & !  7
375                                      15.0_wp,       15.0_wp, 0.00_wp, & !  8
376                                      10.0_wp,       10.0_wp, 0.05_wp, & !  9
377                                      10.0_wp,       10.0_wp, 0.05_wp, & ! 10
378                                      10.0_wp,       10.0_wp, 0.05_wp, & ! 11
379                                      58.0_wp,       58.0_wp, 0.00_wp, & ! 12
380                                      10.0_wp,       10.0_wp, 0.05_wp, & ! 13
381                                    1.0E20_wp,     1.0E20_wp, 0.00_wp, & ! 14
382                                    1.0E20_wp,     1.0E20_wp, 0.00_wp, & ! 15
383                                      10.0_wp,       10.0_wp, 0.05_wp, & ! 16
384                                      10.0_wp,       10.0_wp, 0.05_wp, & ! 17
385                                      20.0_wp,       15.0_wp, 0.03_wp, & ! 18
386                                      20.0_wp,       15.0_wp, 0.03_wp  & ! 19
387                                      /), (/ 3, 19 /) )
388
389!
390!-- Root distribution (sum = 1)  level 1, level 2, level 3, level 4,
391    REAL(wp), DIMENSION(0:3,1:19) :: root_distribution = RESHAPE( (/ &
392                                 0.24_wp, 0.41_wp, 0.31_wp, 0.04_wp, & !  1
393                                 0.35_wp, 0.38_wp, 0.23_wp, 0.04_wp, & !  2
394                                 0.26_wp, 0.39_wp, 0.29_wp, 0.06_wp, & !  3
395                                 0.26_wp, 0.38_wp, 0.29_wp, 0.07_wp, & !  4
396                                 0.24_wp, 0.38_wp, 0.31_wp, 0.07_wp, & !  5
397                                 0.25_wp, 0.34_wp, 0.27_wp, 0.14_wp, & !  6
398                                 0.27_wp, 0.27_wp, 0.27_wp, 0.09_wp, & !  7
399                                 1.00_wp, 0.00_wp, 0.00_wp, 0.00_wp, & !  8
400                                 0.47_wp, 0.45_wp, 0.08_wp, 0.00_wp, & !  9
401                                 0.24_wp, 0.41_wp, 0.31_wp, 0.04_wp, & ! 10
402                                 0.17_wp, 0.31_wp, 0.33_wp, 0.19_wp, & ! 11
403                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 12
404                                 0.25_wp, 0.34_wp, 0.27_wp, 0.11_wp, & ! 13
405                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 14
406                                 0.00_wp, 0.00_wp, 0.00_wp, 0.00_wp, & ! 15
407                                 0.23_wp, 0.36_wp, 0.30_wp, 0.11_wp, & ! 16
408                                 0.23_wp, 0.36_wp, 0.30_wp, 0.11_wp, & ! 17
409                                 0.19_wp, 0.35_wp, 0.36_wp, 0.10_wp, & ! 18
410                                 0.19_wp, 0.35_wp, 0.36_wp, 0.10_wp  & ! 19
411                                 /), (/ 4, 19 /) )
412
413!
414!-- Soil parameters according to the following porosity classes (soil_type)
415
416!
417!-- Soil parameters I           alpha_vg,      l_vg,    n_vg, gamma_w_sat
418    REAL(wp), DIMENSION(0:3,1:7) :: soil_pars = RESHAPE( (/                &
419                                 3.83_wp,  1.250_wp, 1.38_wp,  6.94E-6_wp, & ! 1
420                                 3.14_wp, -2.342_wp, 1.28_wp,  1.16E-6_wp, & ! 2
421                                 0.83_wp, -0.588_wp, 1.25_wp,  0.26E-6_wp, & ! 3
422                                 3.67_wp, -1.977_wp, 1.10_wp,  2.87E-6_wp, & ! 4
423                                 2.65_wp,  2.500_wp, 1.10_wp,  1.74E-6_wp, & ! 5
424                                 1.30_wp,  0.400_wp, 1.20_wp,  0.93E-6_wp, & ! 6
425                                 0.00_wp,  0.00_wp,  0.00_wp,  0.57E-6_wp  & ! 7
426                                 /), (/ 4, 7 /) )
427
428!
429!-- Soil parameters II              m_sat,     m_fc,   m_wilt,    m_res 
430    REAL(wp), DIMENSION(0:3,1:7) :: m_soil_pars = RESHAPE( (/            &
431                                 0.403_wp, 0.244_wp, 0.059_wp, 0.025_wp, & ! 1
432                                 0.439_wp, 0.347_wp, 0.151_wp, 0.010_wp, & ! 2
433                                 0.430_wp, 0.383_wp, 0.133_wp, 0.010_wp, & ! 3
434                                 0.520_wp, 0.448_wp, 0.279_wp, 0.010_wp, & ! 4
435                                 0.614_wp, 0.541_wp, 0.335_wp, 0.010_wp, & ! 5
436                                 0.766_wp, 0.663_wp, 0.267_wp, 0.010_wp, & ! 6
437                                 0.472_wp, 0.323_wp, 0.171_wp, 0.000_wp  & ! 7
438                                 /), (/ 4, 7 /) )
439
440
441    SAVE
442
443
444    PRIVATE
445
446
447!
448!-- Public parameters, constants and initial values
449    PUBLIC alpha_vangenuchten, c_surface, canopy_resistance_coefficient,       &
450           conserve_water_content, dewfall, field_capacity,                    &
451           f_shortwave_incoming, hydraulic_conductivity, init_lsm,             &
452           init_lsm_arrays, lambda_surface_stable, lambda_surface_unstable,    &
453           land_surface, leaf_area_index, lsm_energy_balance, lsm_soil_model,  &
454           lsm_swap_timelevel, l_vangenuchten, min_canopy_resistance,          &
455           min_soil_resistance, n_vangenuchten, residual_moisture, rho_cp,     &
456           rho_lv, root_fraction, saturation_moisture, soil_moisture,          &
457           soil_temperature, soil_type, soil_type_name, vegetation_coverage,   &
458           veg_type, veg_type_name, wilting_point, z0_eb, z0h_eb
459
460!
461!-- Public grid and NetCDF variables
462    PUBLIC dots_soil, id_dim_zs_xy, id_dim_zs_xz, id_dim_zs_yz,                &
463           id_dim_zs_3d, id_dim_zs_mask, id_var_zs_xy, id_var_zs_xz,           &
464           id_var_zs_yz, id_var_zs_3d, id_var_zs_mask, nzb_soil, nzs, nzt_soil,&
465           zs
466
467
468!
469!-- Public 2D output variables
470    PUBLIC c_liq, c_liq_av, c_soil_av, c_veg, c_veg_av, ghf_eb, ghf_eb_av,     &
471           lai, lai_av, qsws_eb, qsws_eb_av, qsws_liq_eb, qsws_liq_eb_av,      &
472           qsws_soil_eb, qsws_soil_eb_av, qsws_veg_eb, qsws_veg_eb_av,         &
473           shf_eb, shf_eb_av
474
475
476!
477!-- Public prognostic variables
478    PUBLIC m_liq_eb, m_liq_eb_av, m_soil, m_soil_av, t_soil, t_soil_av
479
480
481    INTERFACE init_lsm
482       MODULE PROCEDURE init_lsm
483    END INTERFACE init_lsm
484
485    INTERFACE lsm_energy_balance
486       MODULE PROCEDURE lsm_energy_balance
487    END INTERFACE lsm_energy_balance
488
489    INTERFACE lsm_soil_model
490       MODULE PROCEDURE lsm_soil_model
491    END INTERFACE lsm_soil_model
492
493    INTERFACE lsm_swap_timelevel
494       MODULE PROCEDURE lsm_swap_timelevel
495    END INTERFACE lsm_swap_timelevel
496
497 CONTAINS
498
499
500!------------------------------------------------------------------------------!
501! Description:
502! ------------
503!-- Allocate land surface model arrays and define pointers
504!------------------------------------------------------------------------------!
505    SUBROUTINE init_lsm_arrays
506   
507
508       IMPLICIT NONE
509
510!
511!--    Allocate surface and soil temperature / humidity
512#if defined( __nopointer )
513       ALLOCATE ( m_liq_eb(nysg:nyng,nxlg:nxrg) )
514       ALLOCATE ( m_liq_eb_p(nysg:nyng,nxlg:nxrg) )
515       ALLOCATE ( m_soil(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
516       ALLOCATE ( m_soil_p(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
517       ALLOCATE ( t_surface(nysg:nyng,nxlg:nxrg) )
518       ALLOCATE ( t_surface_p(nysg:nyng,nxlg:nxrg) )
519       ALLOCATE ( t_soil(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
520       ALLOCATE ( t_soil_p(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
521#else
522       ALLOCATE ( m_liq_eb_1(nysg:nyng,nxlg:nxrg) )
523       ALLOCATE ( m_liq_eb_2(nysg:nyng,nxlg:nxrg) )
524       ALLOCATE ( m_soil_1(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
525       ALLOCATE ( m_soil_2(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
526       ALLOCATE ( t_surface_1(nysg:nyng,nxlg:nxrg) )
527       ALLOCATE ( t_surface_2(nysg:nyng,nxlg:nxrg) )
528       ALLOCATE ( t_soil_1(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
529       ALLOCATE ( t_soil_2(nzb_soil:nzt_soil+1,nysg:nyng,nxlg:nxrg) )
530#endif
531
532!
533!--    Allocate intermediate timestep arrays
534       ALLOCATE ( tm_liq_eb_m(nysg:nyng,nxlg:nxrg) )
535       ALLOCATE ( tm_soil_m(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
536       ALLOCATE ( tt_surface_m(nysg:nyng,nxlg:nxrg) )
537       ALLOCATE ( tt_soil_m(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
538
539!
540!--    Allocate 2D vegetation model arrays
541       ALLOCATE ( alpha_vg(nysg:nyng,nxlg:nxrg) )
542       ALLOCATE ( c_liq(nysg:nyng,nxlg:nxrg) )
543       ALLOCATE ( c_veg(nysg:nyng,nxlg:nxrg) )
544       ALLOCATE ( f_sw_in(nysg:nyng,nxlg:nxrg) )
545       ALLOCATE ( ghf_eb(nysg:nyng,nxlg:nxrg) )
546       ALLOCATE ( gamma_w_sat(nysg:nyng,nxlg:nxrg) )
547       ALLOCATE ( g_d(nysg:nyng,nxlg:nxrg) )
548       ALLOCATE ( lai(nysg:nyng,nxlg:nxrg) )
549       ALLOCATE ( l_vg(nysg:nyng,nxlg:nxrg) )
550       ALLOCATE ( lambda_h_sat(nysg:nyng,nxlg:nxrg) )
551       ALLOCATE ( lambda_surface_u(nysg:nyng,nxlg:nxrg) )
552       ALLOCATE ( lambda_surface_s(nysg:nyng,nxlg:nxrg) )
553       ALLOCATE ( m_fc(nysg:nyng,nxlg:nxrg) )
554       ALLOCATE ( m_res(nysg:nyng,nxlg:nxrg) )
555       ALLOCATE ( m_sat(nysg:nyng,nxlg:nxrg) )
556       ALLOCATE ( m_wilt(nysg:nyng,nxlg:nxrg) )
557       ALLOCATE ( n_vg(nysg:nyng,nxlg:nxrg) )
558       ALLOCATE ( qsws_eb(nysg:nyng,nxlg:nxrg) )
559       ALLOCATE ( qsws_soil_eb(nysg:nyng,nxlg:nxrg) )
560       ALLOCATE ( qsws_liq_eb(nysg:nyng,nxlg:nxrg) )
561       ALLOCATE ( qsws_veg_eb(nysg:nyng,nxlg:nxrg) )
562       ALLOCATE ( r_a(nysg:nyng,nxlg:nxrg) )
563       ALLOCATE ( r_canopy(nysg:nyng,nxlg:nxrg) )
564       ALLOCATE ( r_soil(nysg:nyng,nxlg:nxrg) )
565       ALLOCATE ( r_soil_min(nysg:nyng,nxlg:nxrg) )
566       ALLOCATE ( r_s(nysg:nyng,nxlg:nxrg) )
567       ALLOCATE ( r_canopy_min(nysg:nyng,nxlg:nxrg) )
568       ALLOCATE ( shf_eb(nysg:nyng,nxlg:nxrg) )
569
570#if ! defined( __nopointer )
571!
572!-- Initial assignment of the pointers
573       t_soil    => t_soil_1;    t_soil_p    => t_soil_2
574       t_surface => t_surface_1; t_surface_p => t_surface_2
575       m_soil    => m_soil_1;    m_soil_p    => m_soil_2
576       m_liq_eb  => m_liq_eb_1;  m_liq_eb_p  => m_liq_eb_2
577#endif
578
579
580    END SUBROUTINE init_lsm_arrays
581
582!------------------------------------------------------------------------------!
583! Description:
584! ------------
585!-- Initialization of the land surface model
586!------------------------------------------------------------------------------!
587    SUBROUTINE init_lsm
588   
589
590       IMPLICIT NONE
591
592       INTEGER(iwp) ::  i !: running index
593       INTEGER(iwp) ::  j !: running index
594       INTEGER(iwp) ::  k !: running index
595
596
597!
598!--    Calculate Exner function
599       exn       = ( surface_pressure / 1000.0_wp )**0.286_wp
600
601
602!
603!--    If no cloud physics is used, rho_surface has not been calculated before
604       IF ( .NOT. cloud_physics )  THEN
605          rho_surface = surface_pressure * 100.0_wp / ( r_d * pt_surface * exn )
606       ENDIF
607
608!
609!--    Calculate frequently used parameters
610       rho_cp    = cp * rho_surface
611       rd_d_rv   = r_d / r_v
612       rho_lv    = rho_surface * l_v
613       drho_l_lv = 1.0_wp / (rho_l * l_v)
614
615!
616!--    Set inital values for prognostic quantities
617       tt_surface_m = 0.0_wp
618       tt_soil_m    = 0.0_wp
619       tm_liq_eb_m  = 0.0_wp
620       c_liq        = 0.0_wp
621
622       ghf_eb = 0.0_wp
623       shf_eb = rho_cp * shf
624
625       IF ( humidity )  THEN
626          qsws_eb = rho_l * l_v * qsws
627       ELSE
628          qsws_eb = 0.0_wp
629       ENDIF
630
631       qsws_liq_eb  = 0.0_wp
632       qsws_soil_eb = qsws_eb
633       qsws_veg_eb  = 0.0_wp
634
635       r_a        = 50.0_wp
636       r_s        = 50.0_wp
637       r_canopy   = 0.0_wp
638       r_soil     = 0.0_wp
639
640!
641!--    Allocate 3D soil model arrays
642       ALLOCATE ( root_fr(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
643       ALLOCATE ( lambda_h(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
644       ALLOCATE ( rho_c_total(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
645
646       lambda_h = 0.0_wp
647!
648!--    If required, allocate humidity-related variables for the soil model
649       IF ( humidity )  THEN
650          ALLOCATE ( lambda_w(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )
651          ALLOCATE ( gamma_w(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) )   
652
653          lambda_w = 0.0_wp 
654       ENDIF
655
656!
657!--    Calculate grid spacings. Temperature and moisture are defined at
658!--    the center of the soil layers, whereas gradients/fluxes are defined
659!--    at the edges (_stag)
660       dz_soil_stag(nzb_soil) = zs(nzb_soil)
661
662       DO k = nzb_soil+1, nzt_soil
663          dz_soil_stag(k) = zs(k) - zs(k-1)
664       ENDDO
665
666       DO k = nzb_soil, nzt_soil-1
667          dz_soil(k) = 0.5 * (dz_soil_stag(k+1) + dz_soil_stag(k))
668       ENDDO
669       dz_soil(nzt_soil) = dz_soil_stag(nzt_soil)
670
671       ddz_soil      = 1.0_wp / dz_soil
672       ddz_soil_stag = 1.0_wp / dz_soil_stag
673
674!
675!--    Initialize standard soil types. It is possible to overwrite each
676!--    parameter by setting the respecticy NAMELIST variable to a
677!--    value /= 9999999.9.
678       IF ( soil_type /= 0 )  THEN 
679 
680          IF ( alpha_vangenuchten == 9999999.9_wp )  THEN
681             alpha_vangenuchten = soil_pars(0,soil_type)
682          ENDIF
683
684          IF ( l_vangenuchten == 9999999.9_wp )  THEN
685             l_vangenuchten = soil_pars(1,soil_type)
686          ENDIF
687
688          IF ( n_vangenuchten == 9999999.9_wp )  THEN
689             n_vangenuchten = soil_pars(2,soil_type)           
690          ENDIF
691
692          IF ( hydraulic_conductivity == 9999999.9_wp )  THEN
693             hydraulic_conductivity = soil_pars(3,soil_type)           
694          ENDIF
695
696          IF ( saturation_moisture == 9999999.9_wp )  THEN
697             saturation_moisture = m_soil_pars(0,soil_type)           
698          ENDIF
699
700          IF ( field_capacity == 9999999.9_wp )  THEN
701             field_capacity = m_soil_pars(1,soil_type)           
702          ENDIF
703
704          IF ( wilting_point == 9999999.9_wp )  THEN
705             wilting_point = m_soil_pars(2,soil_type)           
706          ENDIF
707
708          IF ( residual_moisture == 9999999.9_wp )  THEN
709             residual_moisture = m_soil_pars(3,soil_type)       
710          ENDIF
711
712       ENDIF   
713
714       alpha_vg      = alpha_vangenuchten
715       l_vg          = l_vangenuchten
716       n_vg          = n_vangenuchten 
717       gamma_w_sat   = hydraulic_conductivity
718       m_sat         = saturation_moisture
719       m_fc          = field_capacity
720       m_wilt        = wilting_point
721       m_res         = residual_moisture
722       r_soil_min    = min_soil_resistance
723
724!
725!--    Initial run actions
726       IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
727
728          t_soil    = 0.0_wp
729          m_liq_eb  = 0.0_wp
730          m_soil    = 0.0_wp
731
732!
733!--       Map user settings of T and q for each soil layer
734!--       (make sure that the soil moisture does not drop below the permanent
735!--       wilting point) -> problems with devision by zero)
736          DO k = nzb_soil, nzt_soil
737             t_soil(k,:,:)    = soil_temperature(k)
738             m_soil(k,:,:)    = MAX(soil_moisture(k),m_wilt(:,:))
739             soil_moisture(k) = MAX(soil_moisture(k),wilting_point)
740          ENDDO
741          t_soil(nzt_soil+1,:,:) = soil_temperature(nzt_soil+1)
742
743!
744!--       Calculate surface temperature
745          t_surface = pt_surface * exn
746
747!
748!--       Set artifical values for ts and us so that r_a has its initial value for
749!--       the first time step
750          DO  i = nxlg, nxrg
751             DO  j = nysg, nyng
752                k = nzb_s_inner(j,i)
753
754!
755!--             Assure that r_a cannot be zero at model start
756                IF ( pt(k+1,j,i) == pt(k,j,i) )  pt(k+1,j,i) = pt(k+1,j,i)     &
757                                                 + 1.0E-10_wp
758
759                us(j,i) = 0.1_wp
760                ts(j,i) = (pt(k+1,j,i) - pt(k,j,i)) / r_a(j,i)
761                shf(j,i) = - us(j,i) * ts(j,i)
762             ENDDO
763          ENDDO
764
765!
766!--    Actions for restart runs
767       ELSE
768
769          DO  i = nxlg, nxrg
770             DO  j = nysg, nyng
771                k = nzb_s_inner(j,i)               
772                t_surface(j,i) = pt(k,j,i) * exn
773             ENDDO
774          ENDDO
775
776       ENDIF
777
778!
779!--    Calculate saturation soil heat conductivity
780       lambda_h_sat(:,:) = lambda_h_sm ** (1.0_wp - m_sat(:,:)) *              &
781                           lambda_h_water ** m_sat(:,:)
782
783
784
785
786       DO k = nzb_soil, nzt_soil
787          root_fr(k,:,:) = root_fraction(k)
788       ENDDO
789
790       IF ( veg_type /= 0 )  THEN
791          IF ( min_canopy_resistance == 9999999.9_wp )  THEN
792             min_canopy_resistance = veg_pars(0,veg_type)
793          ENDIF
794          IF ( leaf_area_index == 9999999.9_wp )  THEN
795             leaf_area_index = veg_pars(1,veg_type)         
796          ENDIF
797          IF ( vegetation_coverage == 9999999.9_wp )  THEN
798             vegetation_coverage = veg_pars(2,veg_type)     
799          ENDIF
800          IF ( canopy_resistance_coefficient == 9999999.9_wp )  THEN
801              canopy_resistance_coefficient= veg_pars(3,veg_type)     
802          ENDIF
803          IF ( lambda_surface_stable == 9999999.9_wp )  THEN
804             lambda_surface_stable = surface_pars(0,veg_type)         
805          ENDIF
806          IF ( lambda_surface_unstable == 9999999.9_wp )  THEN
807             lambda_surface_unstable = surface_pars(1,veg_type)       
808          ENDIF
809          IF ( f_shortwave_incoming == 9999999.9_wp )  THEN
810             f_shortwave_incoming = surface_pars(2,veg_type)       
811          ENDIF
812          IF ( z0_eb == 9999999.9_wp )  THEN
813             roughness_length = roughness_par(0,veg_type) 
814             z0_eb            = roughness_par(0,veg_type) 
815          ENDIF
816          IF ( z0h_eb == 9999999.9_wp )  THEN
817             z0h_eb = roughness_par(1,veg_type)
818          ENDIF
819          z0h_factor = z0h_eb / z0_eb
820
821          IF ( ANY( root_fraction == 9999999.9_wp ) )  THEN
822             DO k = nzb_soil, nzt_soil
823                root_fr(k,:,:) = root_distribution(k,veg_type)
824                root_fraction(k) = root_distribution(k,veg_type)
825             ENDDO
826          ENDIF
827
828       ELSE
829
830          IF ( z0_eb == 9999999.9_wp )  THEN
831             z0_eb = roughness_length
832          ENDIF
833          IF ( z0h_eb == 9999999.9_wp )  THEN
834             z0h_eb = z0_eb * z0h_factor
835          ENDIF
836
837       ENDIF
838
839!
840!--    Initialize vegetation
841       r_canopy_min         = min_canopy_resistance
842       lai                  = leaf_area_index
843       c_veg                = vegetation_coverage
844       g_d                  = canopy_resistance_coefficient
845       lambda_surface_s     = lambda_surface_stable
846       lambda_surface_u     = lambda_surface_unstable
847       f_sw_in              = f_shortwave_incoming
848       z0                   = z0_eb
849       z0h                  = z0h_eb
850
851!
852!--    Possibly do user-defined actions (e.g. define heterogeneous land surface)
853       CALL user_init_land_surface
854
855
856       t_soil_p    = t_soil
857       m_soil_p    = m_soil
858       m_liq_eb_p  = m_liq_eb
859
860!--    Store initial profiles of t_soil and m_soil (assuming they are
861!--    horizontally homogeneous on this PE)
862       hom(nzb_soil:nzt_soil,1,90,:)  = SPREAD( t_soil(nzb_soil:nzt_soil,      &
863                                                nysg,nxlg), 2,                 &
864                                                statistic_regions+1 )
865       hom(nzb_soil:nzt_soil,1,92,:)  = SPREAD( m_soil(nzb_soil:nzt_soil,      &
866                                                nysg,nxlg), 2,                 &
867                                                statistic_regions+1 )
868
869!
870!--    Calculate humidity at the surface
871       IF ( humidity )  THEN
872          CALL calc_q_surface
873       ENDIF
874
875
876
877!
878!--    Add timeseries for land surface model
879       dots_label(dots_num+1) = "ghf_eb"
880       dots_label(dots_num+2) = "shf_eb"
881       dots_label(dots_num+3) = "qsws_eb"
882       dots_label(dots_num+4) = "qsws_liq_eb"
883       dots_label(dots_num+5) = "qsws_soil_eb"
884       dots_label(dots_num+6) = "qsws_veg_eb"
885       dots_unit(dots_num+1:dots_num+6) = "W/m2"
886
887       dots_soil = dots_num + 1
888       dots_num  = dots_num + 6
889
890
891       RETURN
892
893    END SUBROUTINE init_lsm
894
895
896
897!------------------------------------------------------------------------------!
898! Description:
899! ------------
900! Solver for the energy balance at the surface.
901!
902! Note: surface fluxes are calculated in the land surface model, but these are
903! not used in the atmospheric code. The fluxes are calculated afterwards in
904! prandtl_fluxes using the surface values of temperature and humidity as
905! provided by the land surface model. In this way, the fluxes in the land
906! surface model are not equal to the ones calculated in prandtl_fluxes
907!------------------------------------------------------------------------------!
908    SUBROUTINE lsm_energy_balance
909
910
911       IMPLICIT NONE
912
913       INTEGER(iwp) ::  i         !: running index
914       INTEGER(iwp) ::  j         !: running index
915       INTEGER(iwp) ::  k, ks     !: running index
916
917       REAL(wp) :: f1,          & !: resistance correction term 1
918                   f2,          & !: resistance correction term 2
919                   f3,          & !: resistance correction term 3
920                   m_min,       & !: minimum soil moisture
921                   T_1,         & !: actual temperature at first grid point
922                   e,           & !: water vapour pressure
923                   e_s,         & !: water vapour saturation pressure
924                   e_s_dT,      & !: derivate of e_s with respect to T
925                   tend,        & !: tendency
926                   dq_s_dT,     & !: derivate of q_s with respect to T
927                   coef_1,      & !: coef. for prognostic equation
928                   coef_2,      & !: coef. for prognostic equation
929                   f_qsws,      & !: factor for qsws_eb
930                   f_qsws_veg,  & !: factor for qsws_veg_eb
931                   f_qsws_soil, & !: factor for qsws_soil_eb
932                   f_qsws_liq,  & !: factor for qsws_liq_eb
933                   f_shf,       & !: factor for shf_eb
934                   lambda_surface, & !: Current value of lambda_surface
935                   m_liq_eb_max   !: maxmimum value of the liq. water reservoir
936
937
938!
939!--    Calculate the exner function for the current time step
940       exn = ( surface_pressure / 1000.0_wp )**0.286_wp
941
942
943       DO i = nxlg, nxrg
944          DO j = nysg, nyng
945             k = nzb_s_inner(j,i)
946
947!
948!--          Set lambda_surface according to stratification
949             IF ( rif(j,i) >= 0.0_wp )  THEN
950                lambda_surface = lambda_surface_s(j,i)
951             ELSE
952                lambda_surface = lambda_surface_u(j,i)
953             ENDIF
954
955!
956!--          First step: calculate aerodyamic resistance. As pt, us, ts
957!--          are not available for the prognostic time step, data from the last
958!--          time step is used here. Note that this formulation is the
959!--          equivalent to the ECMWF formulation using drag coefficients
960             r_a(j,i) = (pt(k+1,j,i) - pt(k,j,i)) / (ts(j,i) * us(j,i) +       &
961                         1.0E-20_wp)
962
963!
964!--          Second step: calculate canopy resistance r_canopy
965!--          f1-f3 here are defined as 1/f1-f3 as in ECMWF documentation
966 
967!--          f1: correction for incoming shortwave radiation (stomata close at
968!--          night)
969             IF ( irad_scheme /= 0 )  THEN
970                f1 = MIN(1.0_wp, ( 0.004_wp * rad_sw_in(j,i) + 0.05_wp ) /     &
971                              (0.81_wp * (0.004_wp * rad_sw_in(j,i) + 1.0_wp) ))
972             ELSE
973                f1 = 1.0_wp
974             ENDIF
975
976!
977!--          f2: correction for soil moisture availability to plants (the
978!--          integrated soil moisture must thus be considered here)
979!--          f2 = 0 for very dry soils
980             m_total = 0.0_wp
981             DO ks = nzb_soil, nzt_soil
982                 m_total = m_total + root_fr(ks,j,i)                           &
983                           * MAX(m_soil(ks,j,i),m_wilt(j,i))
984             ENDDO 
985
986             IF ( m_total .GT. m_wilt(j,i) .AND. m_total .LT. m_fc(j,i) )  THEN
987                f2 = ( m_total - m_wilt(j,i) ) / (m_fc(j,i) - m_wilt(j,i) )
988             ELSEIF ( m_total .GE. m_fc(j,i) )  THEN
989                f2 = 1.0_wp
990             ELSE
991                f2 = 1.0E-20_wp
992             ENDIF
993
994!
995!--          Calculate water vapour pressure at saturation
996             e_s = 0.01_wp * 610.78_wp * EXP( 17.269_wp * ( t_surface(j,i)     &
997                           - 273.16_wp ) / ( t_surface(j,i) - 35.86_wp ) )
998
999!
1000!--          f3: correction for vapour pressure deficit
1001             IF ( g_d(j,i) /= 0.0_wp )  THEN
1002!
1003!--             Calculate vapour pressure
1004                e  = q_p(k+1,j,i) * surface_pressure / 0.622_wp
1005                f3 = EXP ( -g_d(j,i) * (e_s - e) )
1006             ELSE
1007                f3 = 1.0_wp
1008             ENDIF
1009
1010!
1011!--          Calculate canopy resistance. In case that c_veg is 0 (bare soils),
1012!--          this calculation is obsolete, as r_canopy is not used below.
1013!--          To do: check for very dry soil -> r_canopy goes to infinity
1014             r_canopy(j,i) = r_canopy_min(j,i) / (lai(j,i) * f1 * f2 * f3      &
1015                                             + 1.0E-20_wp)
1016
1017!
1018!--          Third step: calculate bare soil resistance r_soil. The Clapp &
1019!--          Hornberger parametrization does not consider c_veg.
1020             IF ( soil_type /= 7 )  THEN
1021                m_min = c_veg(j,i) * m_wilt(j,i) + (1.0_wp - c_veg(j,i)) *     &
1022                        m_res(j,i)
1023             ELSE
1024                m_min = m_wilt(j,i)
1025             ENDIF
1026
1027             f2 = ( m_soil(nzb_soil,j,i) - m_min ) / ( m_fc(j,i) - m_min )
1028             f2 = MAX(f2,1.0E-20_wp)
1029             f2 = MIN(f2,1.0_wp)
1030
1031             r_soil(j,i) = r_soil_min(j,i) / f2
1032
1033!
1034!--          Calculate fraction of liquid water reservoir
1035             m_liq_eb_max = m_max_depth * lai(j,i)
1036             c_liq(j,i) = MIN(1.0_wp, m_liq_eb(j,i) / (m_liq_eb_max+1.0E-20_wp))
1037             
1038
1039!
1040!--          Calculate saturation specific humidity
1041             q_s = 0.622_wp * e_s / surface_pressure
1042
1043!
1044!--          In case of dewfall, set evapotranspiration to zero
1045!--          All super-saturated water is then removed from the air
1046             IF ( humidity .AND. dewfall .AND. q_s .LE. q_p(k+1,j,i) )  THEN
1047                r_canopy(j,i) = 0.0_wp
1048                r_soil(j,i)   = 0.0_wp
1049             ENDIF 
1050
1051!
1052!--          Calculate coefficients for the total evapotranspiration
1053             f_qsws_veg  = rho_lv * c_veg(j,i) * (1.0_wp - c_liq(j,i))/        &
1054                           (r_a(j,i) + r_canopy(j,i))
1055
1056             f_qsws_soil = rho_lv * (1.0_wp - c_veg(j,i)) / (r_a(j,i) +        &
1057                                                          r_soil(j,i))
1058             f_qsws_liq  = rho_lv * c_veg(j,i) * c_liq(j,i) / r_a(j,i)
1059
1060
1061!
1062!--          If soil moisture is below wilting point, plants do no longer
1063!--          transpirate.
1064!              IF ( m_soil(k,j,i) .LT. m_wilt(j,i) )  THEN
1065!                 f_qsws_veg = 0.0_wp
1066!              ENDIF
1067
1068!
1069!--          If dewfall is deactivated, vegetation, soil and liquid water
1070!--          reservoir are not allowed to take up water from the super-saturated
1071!--          air.
1072             IF ( humidity )  THEN
1073                IF ( q_s .LE. q_p(k+1,j,i) )  THEN
1074                   IF ( .NOT. dewfall )  THEN
1075                      f_qsws_veg  = 0.0_wp
1076                      f_qsws_soil = 0.0_wp
1077                      f_qsws_liq  = 0.0_wp
1078                   ENDIF
1079                ENDIF
1080             ENDIF
1081
1082             f_shf  = rho_cp / r_a(j,i)
1083             f_qsws = f_qsws_veg + f_qsws_soil + f_qsws_liq
1084
1085!
1086!--          Calculate derivative of q_s for Taylor series expansion
1087             e_s_dT = e_s * ( 17.269_wp / (t_surface(j,i) - 35.86_wp) -        &
1088                              17.269_wp*(t_surface(j,i) - 273.16_wp)           &
1089                              / (t_surface(j,i) - 35.86_wp)**2 )
1090
1091             dq_s_dT = 0.622_wp * e_s_dT / surface_pressure
1092
1093             T_1 = pt_p(k+1,j,i) * exn
1094
1095!
1096!--          Add LW up so that it can be removed in prognostic equation
1097             rad_net(j,i) = rad_net(j,i) + sigma_SB * t_surface(j,i) ** 4
1098
1099             IF ( humidity )  THEN
1100
1101!
1102!--             Numerator of the prognostic equation
1103                coef_1 = rad_net(j,i) + 3.0_wp * sigma_SB * t_surface(j,i) ** 4&
1104                         + f_shf / exn * T_1 + f_qsws * ( q_p(k+1,j,i) - q_s   &
1105                         + dq_s_dT * t_surface(j,i) ) + lambda_surface         &
1106                         * t_soil(nzb_soil,j,i)
1107
1108!
1109!--             Denominator of the prognostic equation
1110                coef_2 = 4.0_wp * sigma_SB * t_surface(j,i) ** 3 + f_qsws      &
1111                         * dq_s_dT + lambda_surface + f_shf / exn
1112
1113             ELSE
1114
1115!
1116!--             Numerator of the prognostic equation
1117                coef_1 = rad_net(j,i) + 3.0_wp * sigma_SB * t_surface(j,i) ** 4&
1118                         + f_shf / exn * T_1 + lambda_surface                  &
1119                         * t_soil(nzb_soil,j,i)
1120
1121!
1122!--             Denominator of the prognostic equation
1123                coef_2 = 4.0_wp * sigma_SB * t_surface(j,i) ** 3               &
1124                         + lambda_surface + f_shf / exn
1125
1126             ENDIF
1127
1128             tend = 0.0_wp
1129
1130!
1131!--          Implicit solution when the surface layer has no heat capacity,
1132!--          otherwise use RK3 scheme.
1133             t_surface_p(j,i) = ( coef_1 * dt_3d * tsc(2) + c_surface *        &
1134                                t_surface(j,i) ) / ( c_surface + coef_2 * dt_3d&
1135                                * tsc(2) ) 
1136!
1137!--          Add RK3 term
1138             t_surface_p(j,i) = t_surface_p(j,i) + dt_3d * tsc(3)              &
1139                                * tt_surface_m(j,i)
1140
1141!
1142!--          Calculate true tendency
1143             tend = (t_surface_p(j,i) - t_surface(j,i) - dt_3d * tsc(3)        &
1144                    * tt_surface_m(j,i)) / (dt_3d  * tsc(2))
1145!
1146!--          Calculate t_surface tendencies for the next Runge-Kutta step
1147             IF ( timestep_scheme(1:5) == 'runge' )  THEN
1148                IF ( intermediate_timestep_count == 1 )  THEN
1149                   tt_surface_m(j,i) = tend
1150                ELSEIF ( intermediate_timestep_count <                         &
1151                         intermediate_timestep_count_max )  THEN
1152                   tt_surface_m(j,i) = -9.5625_wp * tend + 5.3125_wp           &
1153                                       * tt_surface_m(j,i)
1154                ENDIF
1155             ENDIF
1156
1157             pt_p(k,j,i) = t_surface_p(j,i) / exn
1158!
1159!--          Calculate fluxes
1160             rad_net(j,i)   = rad_net(j,i) + 3.0_wp * sigma_SB                 &
1161                              * t_surface(j,i)**4 - 4.0_wp * sigma_SB          &
1162                              * t_surface(j,i)**3 * t_surface_p(j,i)
1163             ghf_eb(j,i)    = lambda_surface * (t_surface_p(j,i)               &
1164                              - t_soil(nzb_soil,j,i))
1165             shf_eb(j,i)    = - f_shf  * ( pt_p(k+1,j,i) - pt_p(k,j,i) )
1166
1167             IF ( humidity )  THEN
1168                qsws_eb(j,i)  = - f_qsws    * ( q_p(k+1,j,i) - q_s + dq_s_dT   &
1169                                * t_surface(j,i) - dq_s_dT * t_surface_p(j,i) )
1170
1171                qsws_veg_eb(j,i)  = - f_qsws_veg  * ( q_p(k+1,j,i) - q_s       &
1172                                    + dq_s_dT * t_surface(j,i) - dq_s_dT       &
1173                                    * t_surface_p(j,i) )
1174
1175                qsws_soil_eb(j,i) = - f_qsws_soil * ( q_p(k+1,j,i) - q_s       &
1176                                    + dq_s_dT * t_surface(j,i) - dq_s_dT       &
1177                                    * t_surface_p(j,i) )
1178
1179                qsws_liq_eb(j,i)  = - f_qsws_liq  * ( q_p(k+1,j,i) - q_s       &
1180                                    + dq_s_dT * t_surface(j,i) - dq_s_dT       &
1181                                    * t_surface_p(j,i) )
1182             ENDIF
1183
1184!
1185!--          Calculate the true surface resistance
1186             IF ( qsws_eb(j,i) .EQ. 0.0_wp )  THEN
1187                r_s(j,i) = 1.0E10_wp
1188             ELSE
1189                r_s(j,i) = - rho_lv * ( q_p(k+1,j,i) - q_s + dq_s_dT           &
1190                           * t_surface(j,i) - dq_s_dT * t_surface_p(j,i) )     &
1191                           / qsws_eb(j,i) - r_a(j,i)
1192             ENDIF
1193
1194!
1195!--          Calculate change in liquid water reservoir due to dew fall or
1196!--          evaporation of liquid water
1197             IF ( humidity )  THEN
1198!
1199!--             If precipitation is activated, add rain water to qsws_liq_eb.
1200!--             precipitation_rate is given in mm.
1201                IF ( precipitation )  THEN
1202                   qsws_liq_eb(j,i) = qsws_liq_eb(j,i)                         &
1203                                       + precipitation_rate(j,i) * 0.001_wp    &
1204                                       * rho_l * l_v
1205                ENDIF
1206!
1207!--             If the air is saturated, check the reservoir water level
1208                IF ( q_s .LE. q_p(k+1,j,i))  THEN
1209!
1210!--                Check if reservoir is full (avoid values > m_liq_eb_max)
1211!--                In that case, qsws_liq_eb goes to qsws_soil_eb. In this
1212!--                case qsws_veg_eb is zero anyway (because c_liq = 1),       
1213!--                so that tend is zero and no further check is needed
1214                   IF ( m_liq_eb(j,i) .EQ. m_liq_eb_max )  THEN
1215                      qsws_soil_eb(j,i) = qsws_soil_eb(j,i)                    &
1216                                           + qsws_liq_eb(j,i)
1217                      qsws_liq_eb(j,i)  = 0.0_wp
1218                   ENDIF
1219
1220!
1221!--                In case qsws_veg_eb becomes negative (unphysical behavior),
1222!--                let the water enter the liquid water reservoir as dew on the
1223!--                plant
1224                   IF ( qsws_veg_eb(j,i) .LT. 0.0_wp )  THEN
1225                      qsws_liq_eb(j,i) = qsws_liq_eb(j,i) + qsws_veg_eb(j,i)
1226                      qsws_veg_eb(j,i) = 0.0_wp
1227                   ENDIF
1228                ENDIF                   
1229 
1230                tend = - qsws_liq_eb(j,i) * drho_l_lv
1231
1232                m_liq_eb_p(j,i) = m_liq_eb(j,i) + dt_3d * ( tsc(2) * tend      &
1233                                                   + tsc(3) * tm_liq_eb_m(j,i) )
1234
1235!
1236!--             Check if reservoir is overfull -> reduce to maximum
1237!--             (conservation of water is violated here)
1238                m_liq_eb_p(j,i) = MIN(m_liq_eb_p(j,i),m_liq_eb_max)
1239
1240!
1241!--             Check if reservoir is empty (avoid values < 0.0)
1242!--             (conservation of water is violated here)
1243                m_liq_eb_p(j,i) = MAX(m_liq_eb_p(j,i),0.0_wp)
1244
1245
1246!
1247!--             Calculate m_liq_eb tendencies for the next Runge-Kutta step
1248                IF ( timestep_scheme(1:5) == 'runge' )  THEN
1249                   IF ( intermediate_timestep_count == 1 )  THEN
1250                      tm_liq_eb_m(j,i) = tend
1251                   ELSEIF ( intermediate_timestep_count <                      &
1252                            intermediate_timestep_count_max )  THEN
1253                      tm_liq_eb_m(j,i) = -9.5625_wp * tend + 5.3125_wp         &
1254                                      * tm_liq_eb_m(j,i)
1255                   ENDIF
1256                ENDIF
1257
1258             ENDIF
1259
1260          ENDDO
1261       ENDDO
1262
1263    END SUBROUTINE lsm_energy_balance
1264
1265
1266!------------------------------------------------------------------------------!
1267! Description:
1268! ------------
1269!
1270! Soil model as part of the land surface model. The model predicts soil
1271! temperature and water content.
1272!------------------------------------------------------------------------------!
1273    SUBROUTINE lsm_soil_model
1274
1275
1276       IMPLICIT NONE
1277
1278       INTEGER(iwp) ::  i   !: running index
1279       INTEGER(iwp) ::  j   !: running index
1280       INTEGER(iwp) ::  k   !: running index
1281
1282       REAL(wp)     :: h_vg !: Van Genuchten coef. h
1283
1284       REAL(wp), DIMENSION(nzb_soil:nzt_soil) :: gamma_temp,  & !: temp. gamma
1285                                               lambda_temp, & !: temp. lambda
1286                                               tend           !: tendency
1287
1288       DO i = nxlg, nxrg   
1289          DO j = nysg, nyng
1290             DO k = nzb_soil, nzt_soil
1291!
1292!--             Calculate volumetric heat capacity of the soil, taking into
1293!--             account water content
1294                rho_c_total(k,j,i) = (rho_c_soil * (1.0_wp - m_sat(j,i))       &
1295                                     + rho_c_water * m_soil(k,j,i))
1296
1297!
1298!--             Calculate soil heat conductivity at the center of the soil
1299!--             layers
1300                ke = 1.0_wp + LOG10(MAX(0.1_wp,m_soil(k,j,i) / m_sat(j,i)))
1301                lambda_temp(k) = ke * (lambda_h_sat(j,i) + lambda_h_dry) +     &
1302                                 lambda_h_dry
1303
1304             ENDDO
1305
1306!
1307!--          Calculate soil heat conductivity (lambda_h) at the _stag level
1308!--          using linear interpolation
1309             DO k = nzb_soil, nzt_soil-1
1310                 
1311                lambda_h(k,j,i) = lambda_temp(k) +                             &
1312                                  ( lambda_temp(k+1) - lambda_temp(k) )        &
1313                                  * 0.5 * dz_soil_stag(k) * ddz_soil(k+1)
1314
1315             ENDDO
1316             lambda_h(nzt_soil,j,i) = lambda_temp(nzt_soil)
1317
1318!
1319!--          Prognostic equation for soil temperature t_soil
1320             tend(:) = 0.0_wp
1321             tend(0) = (1.0_wp/rho_c_total(nzb_soil,j,i)) *                    &
1322                       ( lambda_h(nzb_soil,j,i) * ( t_soil(nzb_soil+1,j,i)     &
1323                         - t_soil(nzb_soil,j,i) ) * ddz_soil(nzb_soil)         &
1324                         + ghf_eb(j,i) ) * ddz_soil_stag(nzb_soil)
1325
1326             DO  k = 1, nzt_soil
1327                tend(k) = (1.0_wp/rho_c_total(k,j,i))                          &
1328                          * (   lambda_h(k,j,i)                                &
1329                              * ( t_soil(k+1,j,i) - t_soil(k,j,i) )            &
1330                              * ddz_soil(k)                                    &
1331                              - lambda_h(k-1,j,i)                              &
1332                              * ( t_soil(k,j,i) - t_soil(k-1,j,i) )            &
1333                              * ddz_soil(k-1)                                  &
1334                            ) * ddz_soil_stag(k)
1335             ENDDO
1336
1337             t_soil_p(nzb_soil:nzt_soil,j,i) = t_soil(nzb_soil:nzt_soil,j,i)   &
1338                                             + dt_3d * ( tsc(2)                &
1339                                             * tend(:) + tsc(3)                &
1340                                             * tt_soil_m(:,j,i) )   
1341
1342!
1343!--          Calculate t_soil tendencies for the next Runge-Kutta step
1344             IF ( timestep_scheme(1:5) == 'runge' )  THEN
1345                IF ( intermediate_timestep_count == 1 )  THEN
1346                   DO  k = nzb_soil, nzt_soil
1347                      tt_soil_m(k,j,i) = tend(k)
1348                   ENDDO
1349                ELSEIF ( intermediate_timestep_count <                         &
1350                         intermediate_timestep_count_max )  THEN
1351                   DO  k = nzb_soil, nzt_soil
1352                      tt_soil_m(k,j,i) = -9.5625_wp * tend(k) + 5.3125_wp      &
1353                                         * tt_soil_m(k,j,i)
1354                   ENDDO
1355                ENDIF
1356             ENDIF
1357
1358
1359             DO k = nzb_soil, nzt_soil
1360
1361!
1362!--             Calculate soil diffusivity at the center of the soil layers
1363                lambda_temp(k) = (- b_ch * gamma_w_sat(j,i) * psi_sat          &
1364                                  / m_sat(j,i) ) * ( MAX(m_soil(k,j,i),        &
1365                                  m_wilt(j,i)) / m_sat(j,i) )**(b_ch + 2.0_wp)
1366
1367!
1368!--             Parametrization of Van Genuchten
1369                IF ( soil_type /= 7 )  THEN
1370!
1371!--                Calculate the hydraulic conductivity after Van Genuchten
1372!--                (1980)
1373                   h_vg = ( ( (m_res(j,i) - m_sat(j,i)) / ( m_res(j,i) -       &
1374                              MAX(m_soil(k,j,i),m_wilt(j,i)) ) )**(n_vg(j,i)   &
1375                              / (n_vg(j,i)-1.0_wp)) - 1.0_wp                   &
1376                          )**(1.0_wp/n_vg(j,i)) / alpha_vg(j,i)
1377
1378                   gamma_temp(k) = gamma_w_sat(j,i) * ( ( (1.0_wp +            &
1379                                   (alpha_vg(j,i)*h_vg)**n_vg(j,i))**(1.0_wp   &
1380                                   -1.0_wp/n_vg(j,i)) - (alpha_vg(j,i)*h_vg    &
1381                                   )**(n_vg(j,i)-1.0_wp))**2 )                 &
1382                                   / ( (1.0_wp + (alpha_vg(j,i)*h_vg           &
1383                                   )**n_vg(j,i))**((1.0_wp - 1.0_wp/n_vg(j,i)) &
1384                                   *(l_vg(j,i) + 2.0_wp)) )
1385
1386!
1387!--             Parametrization of Clapp & Hornberger
1388                ELSE
1389                   gamma_temp(k) = gamma_w_sat(j,i) * (m_soil(k,j,i)           &
1390                                   / m_sat(j,i) )**(2.0_wp * b_ch + 3.0_wp)
1391                ENDIF
1392
1393             ENDDO
1394
1395
1396             IF ( humidity )  THEN
1397!
1398!--             Calculate soil diffusivity (lambda_w) at the _stag level
1399!--             using linear interpolation. To do: replace this with
1400!--             ECMWF-IFS Eq. 8.81
1401                DO k = nzb_soil, nzt_soil-1
1402                     
1403                   lambda_w(k,j,i) = lambda_temp(k) +                          &
1404                                     ( lambda_temp(k+1) - lambda_temp(k) )     &
1405                                     * 0.5_wp * dz_soil_stag(k) * ddz_soil(k+1)
1406                   gamma_w(k,j,i)  = gamma_temp(k) +                           &
1407                                     ( gamma_temp(k+1) - gamma_temp(k) )       &
1408                                     * 0.5_wp * dz_soil_stag(k) * ddz_soil(k+1)                 
1409
1410                ENDDO
1411
1412!
1413!
1414!--             In case of a closed bottom (= water content is conserved), set
1415!--             hydraulic conductivity to zero to that no water will be lost
1416!--             in the bottom layer.
1417                IF ( conserve_water_content )  THEN
1418                   gamma_w(nzt_soil,j,i) = 0.0_wp
1419                ELSE
1420                   gamma_w(nzt_soil,j,i) = lambda_temp(nzt_soil)
1421                ENDIF     
1422
1423!--             The root extraction (= root_extr * qsws_veg_eb / (rho_l * l_v))
1424!--             ensures the mass conservation for water. The transpiration of
1425!--             plants equals the cumulative withdrawals by the roots in the
1426!--             soil. The scheme takes into account the availability of water
1427!--             in the soil layers as well as the root fraction in the
1428!--             respective layer. Layer with moisture below wilting point will
1429!--             not contribute, which reflects the preference of plants to
1430!--             take water from moister layers.
1431
1432!
1433!--             Calculate the root extraction (ECMWF 7.69, modified so that the
1434!--             sum of root_extr = 1). The energy balance solver guarantees a
1435!--             positive transpiration, so that there is no need for an
1436!--             additional check.
1437
1438                m_total = 0.0_wp
1439                DO k = nzb_soil, nzt_soil
1440                    IF ( m_soil(k,j,i) .GT. m_wilt(j,i) )  THEN
1441                       m_total = m_total + root_fr(k,j,i) * m_soil(k,j,i)
1442                    ENDIF
1443                ENDDO 
1444
1445                DO k = nzb_soil, nzt_soil 
1446                   IF ( m_soil(k,j,i) .GT. m_wilt(j,i) )  THEN
1447                      root_extr(k) = root_fr(k,j,i) * m_soil(k,j,i) / m_total
1448                   ELSE
1449                      root_extr(k) = 0.0_wp
1450                   ENDIF
1451                ENDDO
1452
1453!
1454!--             Prognostic equation for soil water content m_soil
1455                tend(:) = 0.0_wp
1456                tend(nzb_soil) = ( lambda_w(nzb_soil,j,i) * (                  &
1457                            m_soil(nzb_soil+1,j,i) - m_soil(nzb_soil,j,i) )    &
1458                            * ddz_soil(nzb_soil) - gamma_w(nzb_soil,j,i) - (   &
1459                            root_extr(nzb_soil) * qsws_veg_eb(j,i)             &
1460                            + qsws_soil_eb(j,i) ) * drho_l_lv )                &
1461                            * ddz_soil_stag(nzb_soil)
1462
1463                DO  k = nzb_soil+1, nzt_soil-1
1464                   tend(k) = ( lambda_w(k,j,i) * ( m_soil(k+1,j,i)             &
1465                               - m_soil(k,j,i) ) * ddz_soil(k) - gamma_w(k,j,i)&
1466                               - lambda_w(k-1,j,i) * (m_soil(k,j,i) -          &
1467                               m_soil(k-1,j,i)) * ddz_soil(k-1)                &
1468                               + gamma_w(k-1,j,i) - (root_extr(k)              &
1469                               * qsws_veg_eb(j,i) * drho_l_lv)                 &
1470                             ) * ddz_soil_stag(k)
1471
1472                ENDDO
1473                tend(nzt_soil) = ( - gamma_w(nzt_soil,j,i)                     &
1474                                        - lambda_w(nzt_soil-1,j,i)             &
1475                                        * (m_soil(nzt_soil,j,i)                &
1476                                        - m_soil(nzt_soil-1,j,i))              &
1477                                        * ddz_soil(nzt_soil-1)                 &
1478                                        + gamma_w(nzt_soil-1,j,i) - (          &
1479                                          root_extr(nzt_soil)                  &
1480                                        * qsws_veg_eb(j,i) * drho_l_lv  )      &
1481                                      ) * ddz_soil_stag(nzt_soil)             
1482
1483                m_soil_p(nzb_soil:nzt_soil,j,i) = m_soil(nzb_soil:nzt_soil,j,i)&
1484                                                + dt_3d * ( tsc(2) * tend(:)   &
1485                                                + tsc(3) * tm_soil_m(:,j,i) )   
1486
1487!
1488!--             Account for dry soils (find a better solution here!)
1489                m_soil_p(:,j,i) = MAX(m_soil_p(:,j,i),0.0_wp)
1490
1491!
1492!--             Calculate m_soil tendencies for the next Runge-Kutta step
1493                IF ( timestep_scheme(1:5) == 'runge' )  THEN
1494                   IF ( intermediate_timestep_count == 1 )  THEN
1495                      DO  k = nzb_soil, nzt_soil
1496                         tm_soil_m(k,j,i) = tend(k)
1497                      ENDDO
1498                   ELSEIF ( intermediate_timestep_count <                      &
1499                            intermediate_timestep_count_max )  THEN
1500                      DO  k = nzb_soil, nzt_soil
1501                         tm_soil_m(k,j,i) = -9.5625_wp * tend(k) + 5.3125_wp   &
1502                                     * tm_soil_m(k,j,i)
1503                      ENDDO
1504                   ENDIF
1505                ENDIF
1506
1507             ENDIF
1508
1509          ENDDO
1510       ENDDO
1511
1512!
1513!--    Calculate surface specific humidity
1514       IF ( humidity )  THEN
1515          CALL calc_q_surface
1516       ENDIF
1517
1518
1519    END SUBROUTINE lsm_soil_model
1520
1521
1522!------------------------------------------------------------------------------!
1523! Description:
1524! ------------
1525! Calculation of specific humidity of the surface layer (surface)
1526!------------------------------------------------------------------------------!
1527    SUBROUTINE calc_q_surface
1528
1529       IMPLICIT NONE
1530
1531       INTEGER :: i              !: running index
1532       INTEGER :: j              !: running index
1533       INTEGER :: k              !: running index
1534       REAL(wp) :: resistance    !: aerodynamic and soil resistance term
1535
1536       DO i = nxlg, nxrg   
1537          DO j = nysg, nyng
1538             k = nzb_s_inner(j,i)
1539
1540!
1541!--          Calculate water vapour pressure at saturation
1542             e_s = 0.01_wp * 610.78_wp * EXP( 17.269_wp * ( t_surface(j,i)     &
1543                                              - 273.16_wp ) / ( t_surface(j,i) &
1544                                              - 35.86_wp ) )
1545
1546!
1547!--          Calculate specific humidity at saturation
1548             q_s = 0.622_wp * e_s / surface_pressure
1549
1550
1551             resistance = r_a(j,i) / (r_a(j,i) + r_s(j,i))
1552
1553!
1554!--          Calculate specific humidity at surface
1555             q_p(k,j,i) = resistance * q_s + (1.0_wp - resistance)             &
1556                          * q_p(k+1,j,i)
1557
1558          ENDDO
1559       ENDDO
1560
1561    END SUBROUTINE calc_q_surface
1562
1563!------------------------------------------------------------------------------!
1564! Description:
1565! ------------
1566! Swapping of timelevels
1567!------------------------------------------------------------------------------!
1568    SUBROUTINE lsm_swap_timelevel ( mod_count )
1569
1570       IMPLICIT NONE
1571
1572       INTEGER, INTENT(IN) :: mod_count
1573
1574#if defined( __nopointer )
1575
1576       t_surface    = t_surface_p
1577       t_soil       = t_soil_p
1578       IF ( humidity )  THEN
1579          m_soil    = m_soil_p
1580          m_liq_eb  = m_liq_eb_p
1581       ENDIF
1582
1583#else
1584   
1585       SELECT CASE ( mod_count )
1586
1587          CASE ( 0 )
1588
1589             t_surface = t_surface_p
1590             t_soil    = t_soil_p
1591             IF ( humidity )  THEN
1592                m_soil    = m_soil_p
1593                m_liq_eb  = m_liq_eb_p
1594             ENDIF
1595
1596
1597          CASE ( 1 )
1598
1599             t_surface  => t_surface_1; t_surface_p  => t_surface_2
1600             t_soil     => t_soil_1;    t_soil_p     => t_soil_2
1601             IF ( humidity )  THEN
1602                m_soil    => m_soil_1;   m_soil_p    => m_soil_2
1603                m_liq_eb  => m_liq_eb_1; m_liq_eb_p  => m_liq_eb_2
1604             ENDIF
1605
1606       END SELECT
1607#endif
1608
1609    END SUBROUTINE lsm_swap_timelevel
1610
1611
1612 END MODULE land_surface_model_mod
Note: See TracBrowser for help on using the repository browser.