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

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

relicts of rif removed

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