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

Last change on this file since 1783 was 1783, checked in by raasch, 8 years ago

NetCDF routines modularized; new parameter netcdf_deflate; further changes in the pmc

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