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

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

last commit documented

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