SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
mode_flake.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
5 ! File %M% from Library %Q%
6 ! Version %I% from %G% extracted: %H%
7 !------------------------------------------------------------------------------
8 
9 MODULE mode_flake
10 
11 !------------------------------------------------------------------------------
12 !
13 ! Description:
14 !
15 ! The main program unit of the lake model FLake,
16 ! containing most of the FLake procedures.
17 ! Most FLake variables and local parameters are declared.
18 !
19 ! FLake (Fresh-water Lake) is a lake model capable of predicting the surface temperature
20 ! in lakes of various depth on the time scales from a few hours to a year.
21 ! The model is based on a two-layer parametric representation of
22 ! the evolving temperature profile, where the structure of the stratified layer between the
23 ! upper mixed layer and the basin bottom, the lake thermocline,
24 ! is described using the concept of self-similarity of the temperature-depth curve.
25 ! The concept was put forward by Kitaigorodskii and Miropolsky (1970)
26 ! to describe the vertical temperature structure of the oceanic seasonal thermocline.
27 ! It has since been successfully used in geophysical applications.
28 ! The concept of self-similarity of the evolving temperature profile
29 ! is also used to describe the vertical structure of the thermally active upper layer
30 ! of bottom sediments and of the ice and snow cover.
31 !
32 ! The lake model incorporates the heat budget equations
33 ! for the four layers in question, viz., snow, ice, water and bottom sediments,
34 ! developed with due regard for the vertically distributed character
35 ! of solar radiation heating.
36 ! The entrainment equation that incorporates the Zilitinkevich (1975) spin-up term
37 ! is used to compute the depth of a convectively-mixed layer.
38 ! A relaxation-type equation is used
39 ! to compute the wind-mixed layer depth in stable and neutral stratification,
40 ! where a multi-limit formulation for the equilibrium mixed-layer depth
41 ! proposed by Zilitinkevich and Mironov (1996)
42 ! accounts for the effects of the earth's rotation, of the surface buoyancy flux
43 ! and of the static stability in the thermocline.
44 ! The equations for the mixed-layer depth are developed with due regard for
45 ! the volumetric character of the radiation heating.
46 ! Simple thermodynamic arguments are invoked to develop
47 ! the evolution equations for the ice thickness and for the snow thickness.
48 ! The heat flux through the water-bottom sediment interface is computed,
49 ! using a parameterization proposed by Golosov et al. (1998).
50 ! The heat flux trough the air-water interface
51 ! (or through the air-ice or air-snow interface)
52 ! is provided by the driving atmospheric model.
53 !
54 ! Empirical constants and parameters of the lake model
55 ! are estimated, using independent empirical and numerical data.
56 ! They should not be re-evaluated when the model is applied to a particular lake.
57 ! The only lake-specific parameters are the lake depth,
58 ! the optical characteristics of lake water,
59 ! the temperature at the bottom of the thermally active layer
60 ! of bottom sediments and the depth of that layer.
61 !
62 ! A detailed description of the lake model is given in
63 ! Mironov, D. V., 2005:
64 ! Parameterization of Lakes in Numerical Weather Prediction.
65 ! Part 1: Description of a Lake Model.
66 ! Manuscript is available from the author.
67 ! Dmitrii Mironov
68 ! German Weather Service, Kaiserleistr. 29/35, D-63067 Offenbach am Main, Germany.
69 ! dmitrii.mironov@dwd.de
70 !
71 ! Lines embraced with "!_tmp" contain temporary parts of the code.
72 ! Lines embraced/marked with "!_dev" may be replaced
73 ! as improved parameterizations are developed and tested.
74 ! Lines embraced/marked with "!_dm" are DM's comments
75 ! that may be helpful to a user.
76 ! Lines embraced/marked with "!_dbg" are used
77 ! for debugging purposes only.
78 !
79 !
80 ! Current Code Owner: DWD, Dmitrii Mironov
81 ! Phone: +49-69-8062 2705
82 ! Fax: +49-69-8062 3721
83 ! E-mail: dmitrii.mironov@dwd.de
84 !
85 ! History:
86 ! Version Date Name
87 ! ---------- ---------- ----
88 ! 1.00 2005/11/17 Dmitrii Mironov
89 ! Initial release
90 ! !VERSION! !DATE! <Your name>
91 ! <Modification comments>
92 !
93 ! Code Description:
94 ! Language: Fortran 90.
95 ! Software Standards: "European Standards for Writing and
96 ! Documenting Exchangeable Fortran 90 Code".
97 !==============================================================================
98 !
99 ! Declarations:
100 !
101 ! Modules used:
102 
103 !USE modd_data_parameters , ONLY : &
104 ! ireals , &! KIND-type parameter for real variables
105 ! iintegers ! KIND-type parameter for "normal" integer variables
106 
107 !==============================================================================
108 
109 !
110 USE yomhook ,ONLY : lhook, dr_hook
111 USE parkind1 ,ONLY : jprb
112 !
113 USE modd_surfex_omp, ONLY : nblock
114 !
115 IMPLICIT NONE
116 
117 !==============================================================================
118 !
119 ! Declarations
120 !
121 ! The variables declared below
122 ! are accessible to all program units of the MODULE flake.
123 ! Some of them should be USEd by the driving routines that call flake routines.
124 ! These are basically the quantities computed by FLake.
125 ! All variables declared below have a suffix "flk".
126 
127 ! FLake variables of type REAL
128 
129 !RJ: provide default unreasonable value to init for 'ifort -fpic -openmp', to avoid ICE
130 REAL,PARAMETER,PRIVATE :: Z_=-HUGE(0.0)
131 
132 ! Temperatures at the previous time step ("p") and the updated temperatures ("n")
133 REAL :: T_mnw_p_flk=Z_, T_mnw_n_flk=Z_ ! Mean temperature of the water column [K]
134 !$OMP THREADPRIVATE(T_mnw_p_flk, T_mnw_n_flk)
135 REAL :: T_snow_p_flk=Z_, T_snow_n_flk=Z_ ! Temperature at the air-snow interface [K]
136 !$OMP THREADPRIVATE(T_snow_p_flk, T_snow_n_flk)
137 REAL :: T_ice_p_flk=Z_, T_ice_n_flk=Z_ ! Temperature at the snow-ice or air-ice interface [K]
138 !$OMP THREADPRIVATE(T_ice_p_flk, T_ice_n_flk)
139 REAL :: T_wML_p_flk=Z_, T_wML_n_flk=Z_ ! Mixed-layer temperature [K]
140 !$OMP THREADPRIVATE(T_wML_p_flk, T_wML_n_flk)
141 REAL :: T_bot_p_flk=Z_, T_bot_n_flk=Z_ ! Temperature at the water-bottom sediment interface [K]
142 !$OMP THREADPRIVATE(T_bot_p_flk, T_bot_n_flk)
143 REAL :: T_B1_p_flk=Z_, T_B1_n_flk=Z_ ! Temperature at the bottom of the upper layer of the sediments [K]
144 !$OMP THREADPRIVATE(T_B1_p_flk, T_B1_n_flk)
145 
146 ! Thickness of various layers at the previous time step ("p") and the updated values ("n")
147 REAL :: h_snow_p_flk=Z_, h_snow_n_flk=Z_ ! Snow thickness [m]*
148 !$OMP THREADPRIVATE(h_snow_p_flk, h_snow_n_flk)
149 REAL :: h_ice_p_flk=Z_, h_ice_n_flk=Z_ ! Ice thickness [m]
150 !$OMP THREADPRIVATE(h_ice_p_flk, h_ice_n_flk)
151 REAL :: h_ML_p_flk=Z_, h_ML_n_flk=Z_ ! Thickness of the mixed-layer [m]
152 !$OMP THREADPRIVATE(h_ML_p_flk, h_ML_n_flk)
153 REAL :: H_B1_p_flk=Z_, H_B1_n_flk=Z_ ! Thickness of the upper layer of bottom sediments [m]
154 !$OMP THREADPRIVATE(H_B1_p_flk, H_B1_n_flk)
155 
156 ! The shape factor(s) at the previous time step ("p") and the updated value(s) ("n")
157 REAL :: C_T_p_flk=Z_, C_T_n_flk=Z_ ! Shape factor (thermocline)
158 !$OMP THREADPRIVATE(C_T_p_flk, C_T_n_flk)
159 REAL :: C_TT_flk=Z_ ! Dimensionless parameter (thermocline)
160 !$OMP THREADPRIVATE(C_TT_flk)
161 REAL :: C_Q_flk=Z_ ! Shape factor with respect to the heat flux (thermocline)
162 !$OMP THREADPRIVATE(C_Q_flk)
163 REAL :: C_I_flk=Z_ ! Shape factor (ice)
164 !$OMP THREADPRIVATE(C_I_flk)
165 REAL :: C_S_flk=Z_ ! Shape factor (snow)
166 !$OMP THREADPRIVATE(C_S_flk)
167 
168 ! Derivatives of the shape functions
169 REAL :: Phi_T_pr0_flk=Z_ ! d\Phi_T(0)/d\zeta (thermocline)
170 !$OMP THREADPRIVATE(Phi_T_pr0_flk)
171 REAL :: Phi_I_pr0_flk=Z_ ! d\Phi_I(0)/d\zeta_I (ice)
172 !$OMP THREADPRIVATE(Phi_I_pr0_flk)
173 REAL :: Phi_I_pr1_flk=Z_ ! d\Phi_I(1)/d\zeta_I (ice)
174 !$OMP THREADPRIVATE(Phi_I_pr1_flk)
175 REAL :: Phi_S_pr0_flk=Z_ ! d\Phi_S(0)/d\zeta_S (snow)
176 !$OMP THREADPRIVATE(Phi_S_pr0_flk)
177 
178 ! Heat and radiation fluxes
179 REAL :: Q_snow_flk=Z_ ! Heat flux through the air-snow interface [W m^{-2}]
180 !$OMP THREADPRIVATE(Q_snow_flk)
181 REAL :: Q_ice_flk=Z_ ! Heat flux through the snow-ice or air-ice interface [W m^{-2}]
182 !$OMP THREADPRIVATE(Q_ice_flk)
183 REAL :: Q_w_flk=Z_ ! Heat flux through the ice-water or air-water interface [W m^{-2}]
184 !$OMP THREADPRIVATE(Q_w_flk)
185 REAL :: Q_bot_flk=Z_ ! Heat flux through the water-bottom sediment interface [W m^{-2}]
186 !$OMP THREADPRIVATE(Q_bot_flk)
187 REAL :: I_atm_flk=Z_ ! Radiation flux at the lower boundary of the atmosphere [W m^{-2}],
188  ! i.e. the incident radiation flux with no regard for the surface albedo.
189 !$OMP THREADPRIVATE(I_atm_flk)
190 REAL :: I_snow_flk=Z_ ! Radiation flux through the air-snow interface [W m^{-2}]
191 !$OMP THREADPRIVATE(I_snow_flk)
192 REAL :: I_ice_flk=Z_ ! Radiation flux through the snow-ice or air-ice interface [W m^{-2}]
193 !$OMP THREADPRIVATE(I_ice_flk)
194 REAL :: I_w_flk=Z_ ! Radiation flux through the ice-water or air-water interface [W m^{-2}]
195 !$OMP THREADPRIVATE(I_w_flk)
196 REAL :: I_h_flk=Z_ ! Radiation flux through the mixed-layer-thermocline interface [W m^{-2}]
197 !$OMP THREADPRIVATE(I_h_flk)
198 REAL :: I_bot_flk=Z_ ! Radiation flux through the water-bottom sediment interface [W m^{-2}]
199 !$OMP THREADPRIVATE(I_bot_flk)
200 REAL :: I_intm_0_h_flk=Z_ ! Mean radiation flux over the mixed layer [W m^{-1}]
201 !$OMP THREADPRIVATE(I_intm_0_h_flk)
202 REAL :: I_intm_h_D_flk=Z_ ! Mean radiation flux over the thermocline [W m^{-1}]
203 !$OMP THREADPRIVATE(I_intm_h_D_flk)
204 REAL :: Q_star_flk=Z_ ! A generalized heat flux scale [W m^{-2}]
205 !$OMP THREADPRIVATE(Q_star_flk)
206 
207 ! Velocity scales
208 REAL :: u_star_w_flk=Z_ ! Friction velocity in the surface layer of lake water [m s^{-1}]
209 !$OMP THREADPRIVATE(u_star_w_flk)
210 REAL :: w_star_sfc_flk=Z_ ! Convective velocity scale,
211  ! using a generalized heat flux scale [m s^{-1}]
212 !$OMP THREADPRIVATE(w_star_sfc_flk)
213 
214 ! The rate of snow accumulation
215 REAL :: dMsnowdt_flk=Z_ ! The rate of snow accumulation [kg m^{-2} s^{-1}]
216 !$OMP THREADPRIVATE(dMsnowdt_flk)
217 
218 !==============================================================================
219 ! Procedures
220 !==============================================================================
221 
222  CONTAINS
223 
224 !==============================================================================
225 ! The codes of the FLake procedures are stored in separate "*.incf" files
226 ! and are included below.
227 !------------------------------------------------------------------------------
228 
229 !==============================================================================
230 ! For SURFEX needs, separate *.incf files are explicitly expanded
231 !==============================================================================
232 !SURFEX include 'flake_radflux.incf'
233 ! File %M% from Library %Q%
234 ! Version %I% from %G% extracted: %H%
235 !------------------------------------------------------------------------------
236 
237 SUBROUTINE flake_radflux(depth_w,albedo,opticpar_water,opticpar_ice,opticpar_snow)
238 
239 !------------------------------------------------------------------------------
240 !
241 ! Description:
242 !
243 ! Computes the radiation fluxes
244 ! at the snow-ice, ice-water, air-water,
245 ! mixed layer-thermocline and water column-bottom sediment interfaces,
246 ! the mean radiation flux over the mixed layer,
247 ! and the mean radiation flux over the thermocline.
248 !
249 !
250 ! Current Code Owner: DWD, Dmitrii Mironov
251 ! Phone: +49-69-8062 2705
252 ! Fax: +49-69-8062 3721
253 ! E-mail: dmitrii.mironov@dwd.de
254 !
255 ! History:
256 ! Version Date Name
257 ! ---------- ---------- ----
258 ! 1.00 2005/11/17 Dmitrii Mironov
259 ! Initial release
260 ! !VERSION! !DATE! <Your name>
261 ! <Modification comments>
262 !
263 ! Code Description:
264 ! Language: Fortran 90.
265 ! Software Standards: "European Standards for Writing and
266 ! Documenting Exchangeable Fortran 90 Code".
267 !==============================================================================
268 !
269 ! Declarations:
270 !
271 ! Modules used:
272 
273 !_dm Parameters are USEd in module "flake".
274 !_nu USE modd_data_parameters , ONLY : &
275 !_nu ireals, & ! KIND-type parameter for real variables
276 !_nu iintegers ! KIND-type parameter for "normal" integer variables
277 
278 USE modd_flake_derivedtypes ! Definitions of derived TYPEs
279 
280 USE modd_flake_parameters , ONLY : &
281  h_snow_min_flk , &! Minimum snow thickness [m]
282  h_ice_min_flk , &! Minimum ice thickness [m]
283  h_ml_min_flk ! Minimum mixed-layer depth [m]
284 
285 !==============================================================================
286 
287 IMPLICIT NONE
288 
289 !==============================================================================
290 !
291 ! Declarations
292 
293 ! Input (procedure arguments)
294 
295 REAL, INTENT(IN) :: &
296  depth_w , &! The lake depth [m]
297  albedo ! Albedo of all surfaces
298 
299 TYPE (opticpar_medium), INTENT(IN) :: &
300  opticpar_water , &! Optical characteristics of water
301  opticpar_ice , &! Optical characteristics of ice
302  opticpar_snow ! Optical characteristics of snow
303 
304 
305 ! Local variables of type INTEGER
306 INTEGER :: &! Help variable(s)
307  i ! DO loop index
308 REAL(KIND=JPRB) :: zhook_handle
309 
310 !==============================================================================
311 ! Start calculations
312 !------------------------------------------------------------------------------
313 
314  IF (lhook) CALL dr_hook('FLAKE:FLAKE_RADFLUX',0,zhook_handle)
315  IF(h_ice_p_flk.GE.h_ice_min_flk) THEN ! Ice exists
316  IF(h_snow_p_flk.GE.h_snow_min_flk) THEN ! There is snow above the ice
317  i_snow_flk = i_atm_flk*(1.-albedo)
318  i_bot_flk = 0.
319  DO i=1, opticpar_snow%nband_optic
320  i_bot_flk = i_bot_flk + &
321  opticpar_snow%frac_optic(i)*exp(-opticpar_snow%extincoef_optic(i)*h_snow_p_flk)
322  END DO
323  i_ice_flk = i_snow_flk*i_bot_flk
324  ELSE ! No snow above the ice
325  i_snow_flk = i_atm_flk
326  i_ice_flk = i_atm_flk*(1.-albedo)
327  END IF
328  i_bot_flk = 0.
329  DO i=1, opticpar_ice%nband_optic
330  i_bot_flk = i_bot_flk + &
331  opticpar_ice%frac_optic(i)*exp(-opticpar_ice%extincoef_optic(i)*h_ice_p_flk)
332  END DO
333  i_w_flk = i_ice_flk*i_bot_flk
334  ELSE ! No ice-snow cover
335  i_snow_flk = i_atm_flk
336  i_ice_flk = i_atm_flk
337  i_w_flk = i_atm_flk*(1.-albedo)
338  END IF
339 
340  IF(h_ml_p_flk.GE.h_ml_min_flk) THEN ! Radiation flux at the bottom of the mixed layer
341  i_bot_flk = 0.
342  DO i=1, opticpar_water%nband_optic
343  i_bot_flk = i_bot_flk + &
344  opticpar_water%frac_optic(i)*exp(-opticpar_water%extincoef_optic(i)*h_ml_p_flk)
345  END DO
346  i_h_flk = i_w_flk*i_bot_flk
347  ELSE ! Mixed-layer depth is less then a minimum value
348  i_h_flk = i_w_flk
349  END IF
350 
351  i_bot_flk = 0. ! Radiation flux at the lake bottom
352  DO i=1, opticpar_water%nband_optic
353  i_bot_flk = i_bot_flk + &
354  opticpar_water%frac_optic(i)*exp(-opticpar_water%extincoef_optic(i)*depth_w)
355  END DO
356  i_bot_flk = i_w_flk*i_bot_flk
357 
358  IF(h_ml_p_flk.GE.h_ml_min_flk) THEN ! Integral-mean radiation flux over the mixed layer
359  i_intm_0_h_flk = 0.
360  DO i=1, opticpar_water%nband_optic
361  i_intm_0_h_flk = i_intm_0_h_flk + &
362  opticpar_water%frac_optic(i)/opticpar_water%extincoef_optic(i)* &
363  (1. - exp(-opticpar_water%extincoef_optic(i)*h_ml_p_flk))
364  END DO
365  i_intm_0_h_flk = i_w_flk*i_intm_0_h_flk/h_ml_p_flk
366  ELSE
367  i_intm_0_h_flk = i_h_flk
368  END IF
369 
370  IF(h_ml_p_flk.LE.depth_w-h_ml_min_flk) THEN ! Integral-mean radiation flux over the thermocline
371  i_intm_h_d_flk = 0.
372  DO i=1, opticpar_water%nband_optic
373  i_intm_h_d_flk = i_intm_h_d_flk + &
374  opticpar_water%frac_optic(i)/opticpar_water%extincoef_optic(i)* &
375  ( exp(-opticpar_water%extincoef_optic(i)*h_ml_p_flk) &
376  - exp(-opticpar_water%extincoef_optic(i)*depth_w) )
377  END DO
378  i_intm_h_d_flk = i_w_flk*i_intm_h_d_flk/(depth_w-h_ml_p_flk)
379  ELSE
380  i_intm_h_d_flk = i_h_flk
381  END IF
382 IF (lhook) CALL dr_hook('FLAKE:FLAKE_RADFLUX',1,zhook_handle)
383 
384 !------------------------------------------------------------------------------
385 ! End calculations
386 !==============================================================================
387 
388 END SUBROUTINE flake_radflux
389 
390 
391 !==============================================================================
392 
393 !==============================================================================
394 !SURFEX include 'flake_driver.incf'
395 ! File %M% from Library %Q%
396 ! Version %I% from %G% extracted: %H%
397 !------------------------------------------------------------------------------
398 
399 SUBROUTINE flake_driver ( depth_w, depth_bs, T_bs, par_Coriolis, &
400  extincoef_water_typ, &
401  del_time, t_sfc_p, t_sfc_n )
402 
403 !------------------------------------------------------------------------------
404 !
405 ! Description:
406 !
407 ! The main driving routine of the lake model FLake
408 ! where computations are performed.
409 ! Advances the surface temperature
410 ! and other FLake variables one time step.
411 ! At the moment, the Euler explicit scheme is used.
412 !
413 ! Lines embraced with "!_tmp" contain temporary parts of the code.
414 ! Lines embraced/marked with "!_dev" may be replaced
415 ! as improved parameterizations are developed and tested.
416 ! Lines embraced/marked with "!_dm" are DM's comments
417 ! that may be helpful to a user.
418 ! Lines embraced/marked with "!_dbg" are used
419 ! for debugging purposes only.
420 !
421 !
422 ! Current Code Owner: DWD, Dmitrii Mironov
423 ! Phone: +49-69-8062 2705
424 ! Fax: +49-69-8062 3721
425 ! E-mail: dmitrii.mironov@dwd.de
426 !
427 ! History:
428 ! Version Date Name
429 ! ---------- ---------- ----
430 ! 1.00 2005/11/17 Dmitrii Mironov
431 ! Initial release
432 ! !VERSION! !DATE! <Your name>
433 ! <Modification comments>
434 !
435 ! Code Description:
436 ! Language: Fortran 90.
437 ! Software Standards: "European Standards for Writing and
438 ! Documenting Exchangeable Fortran 90 Code".
439 !==============================================================================
440 !
441 ! Declarations:
442 !
443 ! Modules used:
444 
445 !_dm Parameters are USEd in module "flake".
446 !_nu USE modd_data_parameters , ONLY : &
447 !_nu ireals, & ! KIND-type parameter for real variables
448 !_nu iintegers ! KIND-type parameter for "normal" integer variables
449 
450 USE modd_flake_parameters ! Thermodynamic parameters and dimensionless constants of FLake
451 
452 USE modd_flake_configure ! Switches and parameters that configure FLake
453 
454 !==============================================================================
455 
456 IMPLICIT NONE
457 
458 !==============================================================================
459 !
460 ! Declarations
461 
462 ! Input (procedure arguments)
463 
464 REAL, INTENT(IN) :: &
465  depth_w , &! The lake depth [m]
466  depth_bs , &! Depth of the thermally active layer of bottom sediments [m]
467  t_bs , &! Temperature at the outer edge of
468  ! the thermally active layer of bottom sediments [K]
469  par_coriolis , &! The Coriolis parameter [s^{-1}]
470  extincoef_water_typ , &! "Typical" extinction coefficient of the lake water [m^{-1}],
471  ! used to compute the equilibrium CBL depth
472  del_time , &! The model time step [s]
473  t_sfc_p ! Surface temperature at the previous time step [K]
474  ! (equal to either T_ice, T_snow or to T_wML)
475 
476 ! Output (procedure arguments)
477 
478 REAL, INTENT(OUT) :: &
479  t_sfc_n ! Updated surface temperature [K]
480  ! (equal to the updated value of either T_ice, T_snow or T_wML)
481 
482 
483 ! Local variables of type LOGICAL
484 LOGICAL :: &
485  l_ice_create , &! Switch, .TRUE. = ice does not exist but should be created
486  l_snow_exists , &! Switch, .TRUE. = there is snow above the ice
487  l_ice_meltabove ! Switch, .TRUE. = snow/ice melting from above takes place
488 
489 ! Local variables of type INTEGER
490 INTEGER :: &
491  i ! Loop index
492 
493 ! Local variables of type REAL
494 REAL :: &
495  d_t_mnw_dt , &! Time derivative of T_mnw [K s^{-1}]
496  d_t_ice_dt , &! Time derivative of T_ice [K s^{-1}]
497  d_t_bot_dt , &! Time derivative of T_bot [K s^{-1}]
498  d_t_b1_dt , &! Time derivative of T_B1 [K s^{-1}]
499  d_h_snow_dt , &! Time derivative of h_snow [m s^{-1}]
500  d_h_ice_dt , &! Time derivative of h_ice [m s^{-1}]
501  d_h_ml_dt , &! Time derivative of h_ML [m s^{-1}]
502  d_h_b1_dt , &! Time derivative of H_B1 [m s^{-1}]
503  d_c_t_dt ! Time derivative of C_T [s^{-1}]
504 
505 ! Local variables of type REAL
506 REAL :: &
507  n_t_mean , &! The mean buoyancy frequency in the thermocline [s^{-1}]
508  zm_h_scale , &! The ZM96 equilibrium SBL depth scale [m]
509  conv_equil_h_scale ! The equilibrium CBL depth scale [m]
510 
511 ! Local variables of type REAL
512 REAL :: &
513  h_ice_threshold , &! If h_ice<h_ice_threshold, use quasi-equilibrium ice model
514  flk_str_1 , &! Help storage variable
515  flk_str_2 , &! Help storage variable
516  r_h_icesnow , &! Dimensionless ratio, used to store intermediate results
517  r_rho_c_icesnow , &! Dimensionless ratio, used to store intermediate results
518  r_ti_icesnow , &! Dimensionless ratio, used to store intermediate results
519  r_tstar_icesnow ! Dimensionless ratio, used to store intermediate results
520 REAL(KIND=JPRB) :: zhook_handle
521 
522 !rsal
523 !REAL :: aux1,aux2,aux3,aux4,aux5
524 
525 !==============================================================================
526 ! Start calculations
527 !------------------------------------------------------------------------------
528 
529 !_dm
530 ! Security. Set time-rate-of-change of prognostic variables to zero.
531 ! Set prognostic variables to their values at the previous time step.
532 ! (This is to avoid spurious changes of prognostic variables
533 ! when FLake is used within a 3D model, e.g. to avoid spurious generation of ice
534 ! at the neighbouring lake points as noticed by Burkhardt Rockel.)
535 !_dm
536 
537 
538 IF (lhook) CALL dr_hook('FLAKE:FLAKE_DRIVER',0,zhook_handle)
539 d_t_mnw_dt = 0.
540 d_t_ice_dt = 0.
541 d_t_bot_dt = 0.
542 d_t_b1_dt = 0.
543 d_h_snow_dt = 0.
544 d_h_ice_dt = 0.
545 d_h_ml_dt = 0.
546 d_h_b1_dt = 0.
547 d_c_t_dt = 0.
548 t_snow_n_flk = t_snow_p_flk
549 t_ice_n_flk = t_ice_p_flk
550 t_wml_n_flk = t_wml_p_flk
551 t_mnw_n_flk = t_mnw_p_flk
552 t_bot_n_flk = t_bot_p_flk
553 t_b1_n_flk = t_b1_p_flk
554 h_snow_n_flk = h_snow_p_flk
555 h_ice_n_flk = h_ice_p_flk
556 h_ml_n_flk = h_ml_p_flk
557 h_b1_n_flk = h_b1_p_flk
558  c_t_n_flk = c_t_p_flk
559 
560 !------------------------------------------------------------------------------
561 ! Compute fluxes, using variables from the previous time step.
562 !------------------------------------------------------------------------------
563 
564 !_dm
565 ! At this point, the heat and radiation fluxes, namely,
566 ! Q_snow_flk, Q_ice_flk, Q_w_flk,
567 ! I_atm_flk, I_snow_flk, I_ice_flk, I_w_flk, I_h_flk, I_bot_flk,
568 ! the mean radiation flux over the mixed layer, I_intm_0_h_flk,
569 ! and the mean radiation flux over the thermocline, I_intm_h_D_flk,
570 ! should be known.
571 ! They are computed within "flake_interface" (or within the driving model)
572 ! and are available to "flake_driver"
573 ! through the above variables declared in the MODULE "flake".
574 ! In case a lake is ice-covered, Q_w_flk is re-computed below.
575 !_dm
576 
577 ! Heat flux through the ice-water interface
578 IF(h_ice_p_flk.GE.h_ice_min_flk) THEN ! Ice exists
579  IF(h_ml_p_flk.LE.h_ml_min_flk) THEN ! Mixed-layer depth is zero, compute flux
580  q_w_flk = -tpl_kappa_w*(t_bot_p_flk-t_wml_p_flk)/depth_w ! Flux with linear T(z)
581  phi_t_pr0_flk = phi_t_pr0_1*c_t_p_flk-phi_t_pr0_2 ! d\Phi(0)/d\zeta (thermocline)
582  q_w_flk = q_w_flk*max(phi_t_pr0_flk, 1.) ! Account for an increased d\Phi(0)/d\zeta
583  ELSE
584  q_w_flk = 0. ! Mixed-layer depth is greater than zero, set flux to zero
585  END IF
586 END IF
587 
588 ! A generalized heat flux scale
589 q_star_flk = q_w_flk + i_w_flk + i_h_flk - 2.*i_intm_0_h_flk
590 
591 ! Heat flux through the water-bottom sediment interface
592 IF(lflk_botsed_use) THEN
593  q_bot_flk = -tpl_kappa_w*(t_b1_p_flk-t_bot_p_flk)/max(h_b1_p_flk, h_b1_min_flk)*phi_b1_pr0
594 ELSE
595  q_bot_flk = 0. ! The bottom-sediment scheme is not used
596 END IF
597 
598 
599 !------------------------------------------------------------------------------
600 ! Check if ice exists or should be created.
601 ! If so, compute the thickness and the temperature of ice and snow.
602 !------------------------------------------------------------------------------
603 
604 !_dm
605 ! Notice that a quasi-equilibrium ice-snow model is used
606 ! to avoid numerical instability when the ice is thin.
607 ! This is always the case when new ice is created.
608 !_dm
609 
610 !_dev
611 ! The dependence of snow density and of snow heat conductivity
612 ! on the snow thickness is accounted for parametrically.
613 ! That is, the time derivatives of \rho_S and \kappa_S are neglected.
614 ! The exception is the equation for the snow thickness
615 ! in case of snow accumulation and no melting,
616 ! where d\rho_S/dt is incorporated.
617 ! Furthermore, some (presumably small) correction terms incorporating
618 ! the snow density and the snow heat conductivity are dropped out.
619 ! Those terms may be included as better formulations
620 ! for \rho_S and \kappa_S are available.
621 !_dev
622 
623 ! Default values
624 l_ice_create = .false.
625 l_ice_meltabove = .false.
626 
627 ice_exist: IF(h_ice_p_flk.LT.h_ice_min_flk) THEN ! Ice does not exist
628 
629  l_ice_create = t_wml_p_flk.LE.(tpl_t_f+c_small_flk).AND.q_w_flk.LT.0.
630  IF(l_ice_create) THEN ! Ice does not exist but should be created
631  d_h_ice_dt = -q_w_flk/tpl_rho_i/tpl_l_f
632  h_ice_n_flk = h_ice_p_flk + d_h_ice_dt*del_time ! Advance h_ice
633  t_ice_n_flk = tpl_t_f + h_ice_n_flk*q_w_flk/tpl_kappa_i/phi_i_pr0_lin ! Ice temperature
634  d_h_snow_dt = dmsnowdt_flk/tpl_rho_s_min
635  h_snow_n_flk = h_snow_p_flk + d_h_snow_dt*del_time ! Advance h_snow
636  phi_i_pr1_flk = phi_i_pr1_lin &
637  + phi_i_ast_mr*min(1., h_ice_n_flk/h_ice_max) ! d\Phi_I(1)/d\zeta_I (ice)
638  r_h_icesnow = phi_i_pr1_flk/phi_s_pr0_lin*tpl_kappa_i/flake_snowheatconduct(h_snow_n_flk) &
639  * h_snow_n_flk/max(h_ice_n_flk, h_ice_min_flk)
640  t_snow_n_flk = t_ice_n_flk + r_h_icesnow*(t_ice_n_flk-tpl_t_f) ! Snow temperature
641  ELSE !kk, 30092004 - to avoid unsertainty here
642  h_ice_n_flk = h_ice_p_flk
643  h_snow_n_flk = h_snow_p_flk
644  t_ice_n_flk = t_ice_p_flk
645  t_snow_n_flk = t_snow_p_flk
646  END IF
647 
648 ELSE ice_exist ! Ice exists
649 
650  l_snow_exists = h_snow_p_flk.GE.h_snow_min_flk ! Check if there is snow above the ice
651 
652  melting: IF(t_snow_p_flk.GE.(tpl_t_f-c_small_flk)) THEN ! T_sfc = T_f, check for melting from above
653  ! T_snow = T_ice if snow is absent
654  IF(l_snow_exists) THEN ! There is snow above the ice
655  flk_str_1 = q_snow_flk + i_snow_flk - i_ice_flk ! Atmospheric forcing
656  IF(flk_str_1.GE.0.) THEN ! Melting of snow and ice from above
657  l_ice_meltabove = .true.
658  d_h_snow_dt = (-flk_str_1/tpl_l_f+dmsnowdt_flk)/flake_snowdensity(h_snow_p_flk)
659  d_h_ice_dt = -(i_ice_flk - i_w_flk - q_w_flk)/tpl_l_f/tpl_rho_i
660  END IF
661  ELSE ! No snow above the ice
662  flk_str_1 = q_ice_flk + i_ice_flk - i_w_flk - q_w_flk ! Atmospheric forcing + heating from the water
663  IF(flk_str_1.GE.0.) THEN ! Melting of ice from above, snow accumulation may occur
664  l_ice_meltabove = .true.
665  d_h_ice_dt = -flk_str_1/tpl_l_f/tpl_rho_i
666  d_h_snow_dt = dmsnowdt_flk/tpl_rho_s_min
667  END IF
668  END IF
669  IF(l_ice_meltabove) THEN ! Melting from above takes place
670  h_ice_n_flk = h_ice_p_flk + d_h_ice_dt *del_time ! Advance h_ice
671  h_snow_n_flk = h_snow_p_flk + d_h_snow_dt*del_time ! Advance h_snow
672  t_ice_n_flk = tpl_t_f ! Set T_ice to the freezing point
673  t_snow_n_flk = tpl_t_f ! Set T_snow to the freezing point
674  END IF
675 
676  END IF melting
677 
678  no_melting: IF(.NOT.l_ice_meltabove) THEN ! No melting from above
679 
680  d_h_snow_dt = flake_snowdensity(h_snow_p_flk)
681  IF(d_h_snow_dt.LT.tpl_rho_s_max) THEN ! Account for d\rho_S/dt
682  flk_str_1 = h_snow_p_flk*tpl_gamma_rho_s/tpl_rho_w_r
683  flk_str_1 = flk_str_1/(1.-flk_str_1)
684  ELSE ! Snow density is equal to its maximum value, d\rho_S/dt=0
685  flk_str_1 = 0.
686  END IF
687  d_h_snow_dt = dmsnowdt_flk/d_h_snow_dt/(1.+flk_str_1) ! Snow accumulation
688  h_snow_n_flk = h_snow_p_flk + d_h_snow_dt*del_time ! Advance h_snow
689 
690  phi_i_pr0_flk = h_ice_p_flk/h_ice_max ! h_ice relative to its maximum value
691  c_i_flk = c_i_lin - c_i_mr*(1.+phi_i_ast_mr)*phi_i_pr0_flk ! Shape factor (ice)
692  phi_i_pr1_flk = phi_i_pr1_lin + phi_i_ast_mr*phi_i_pr0_flk ! d\Phi_I(1)/d\zeta_I (ice)
693  phi_i_pr0_flk = phi_i_pr0_lin - phi_i_pr0_flk ! d\Phi_I(0)/d\zeta_I (ice)
694 
695  h_ice_threshold = max(1., 2.*c_i_flk*tpl_c_i*(tpl_t_f-t_ice_p_flk)/tpl_l_f)
696  h_ice_threshold = phi_i_pr0_flk/c_i_flk*tpl_kappa_i/tpl_rho_i/tpl_c_i*h_ice_threshold
697  h_ice_threshold = sqrt(h_ice_threshold*del_time) ! Threshold value of h_ice
698  h_ice_threshold = min(0.9*h_ice_max, max(h_ice_threshold, h_ice_min_flk))
699  ! h_ice(threshold) < 0.9*H_Ice_max
700 
701  IF(h_ice_p_flk.LT.h_ice_threshold) THEN ! Use a quasi-equilibrium ice model
702 
703  IF(l_snow_exists) THEN ! Use fluxes at the air-snow interface
704  flk_str_1 = q_snow_flk + i_snow_flk - i_w_flk
705  ELSE ! Use fluxes at the air-ice interface
706  flk_str_1 = q_ice_flk + i_ice_flk - i_w_flk
707  END IF
708  d_h_ice_dt = -(flk_str_1-q_w_flk)/tpl_l_f/tpl_rho_i
709  h_ice_n_flk = h_ice_p_flk + d_h_ice_dt *del_time ! Advance h_ice
710  t_ice_n_flk = tpl_t_f + h_ice_n_flk*flk_str_1/tpl_kappa_i/phi_i_pr0_flk ! Ice temperature
711 
712  ELSE ! Use a complete ice model
713 
714  d_h_ice_dt = tpl_kappa_i*(tpl_t_f-t_ice_p_flk)/h_ice_p_flk*phi_i_pr0_flk
715  d_h_ice_dt = (q_w_flk+d_h_ice_dt)/tpl_l_f/tpl_rho_i
716  h_ice_n_flk = h_ice_p_flk + d_h_ice_dt*del_time ! Advance h_ice
717 
718  r_ti_icesnow = tpl_c_i*(tpl_t_f-t_ice_p_flk)/tpl_l_f ! Dimensionless parameter
719  r_tstar_icesnow = 1. - c_i_flk ! Dimensionless parameter
720  IF(l_snow_exists) THEN ! There is snow above the ice
721  r_h_icesnow = phi_i_pr1_flk/phi_s_pr0_lin*tpl_kappa_i/flake_snowheatconduct(h_snow_p_flk) &
722  * h_snow_p_flk/h_ice_p_flk
723  r_rho_c_icesnow = flake_snowdensity(h_snow_p_flk)*tpl_c_s/tpl_rho_i/tpl_c_i
724 !_dev
725 !_dm
726 ! These terms should be included as an improved understanding of the snow scheme is gained,
727 ! of the effect of snow density in particular.
728 !_dm
729 !_nu R_Tstar_icesnow = R_Tstar_icesnow &
730 !_nu + (1.+C_S_lin*h_snow_p_flk/h_ice_p_flk)*R_H_icesnow*R_rho_c_icesnow
731 !_dev
732 
733  r_tstar_icesnow = r_tstar_icesnow*r_ti_icesnow ! Dimensionless parameter
734 
735 !_dev
736 !_nu R_Tstar_icesnow = R_Tstar_icesnow &
737 !_nu + (1.-R_rho_c_icesnow)*tpl_c_I*T_ice_p_flk/tpl_L_f
738 !_dev
739  flk_str_2 = q_snow_flk+i_snow_flk-i_w_flk ! Atmospheric fluxes
740  flk_str_1 = c_i_flk*h_ice_p_flk + (1.+c_s_lin*r_h_icesnow)*r_rho_c_icesnow*h_snow_p_flk
741  d_t_ice_dt = -(1.-2.*c_s_lin)*r_h_icesnow*(tpl_t_f-t_ice_p_flk) &
742  * tpl_c_s*dmsnowdt_flk ! Effect of snow accumulation
743  ELSE ! No snow above the ice
744  r_tstar_icesnow = r_tstar_icesnow*r_ti_icesnow ! Dimensionless parameter
745  flk_str_2 = q_ice_flk+i_ice_flk-i_w_flk ! Atmospheric fluxes
746  flk_str_1 = c_i_flk*h_ice_p_flk
747  d_t_ice_dt = 0.
748  END IF
749  d_t_ice_dt = d_t_ice_dt + tpl_kappa_i*(tpl_t_f-t_ice_p_flk)/h_ice_p_flk*phi_i_pr0_flk &
750  * (1.-r_tstar_icesnow) ! Add flux due to heat conduction
751  d_t_ice_dt = d_t_ice_dt - r_tstar_icesnow*q_w_flk ! Add flux from water to ice
752  d_t_ice_dt = d_t_ice_dt + flk_str_2 ! Add atmospheric fluxes
753  d_t_ice_dt = d_t_ice_dt/tpl_rho_i/tpl_c_i ! Total forcing
754  d_t_ice_dt = d_t_ice_dt/flk_str_1 ! dT_ice/dt
755  t_ice_n_flk = t_ice_p_flk + d_t_ice_dt*del_time ! Advance T_ice
756  END IF
757 
758  phi_i_pr1_flk = min(1., h_ice_n_flk/h_ice_max) ! h_ice relative to its maximum value
759  phi_i_pr1_flk = phi_i_pr1_lin + phi_i_ast_mr*phi_i_pr1_flk ! d\Phi_I(1)/d\zeta_I (ice)
760  r_h_icesnow = phi_i_pr1_flk/phi_s_pr0_lin*tpl_kappa_i/flake_snowheatconduct(h_snow_n_flk) &
761  *h_snow_n_flk/max(h_ice_n_flk, h_ice_min_flk)
762  t_snow_n_flk = t_ice_n_flk + r_h_icesnow*(t_ice_n_flk-tpl_t_f) ! Snow temperature
763 
764  END IF no_melting
765 
766 END IF ice_exist
767 
768 ! Security, limit h_ice by its maximum value
769 h_ice_n_flk = min(h_ice_n_flk, h_ice_max)
770 
771 ! Security, limit the ice and snow temperatures by the freezing point
772 t_snow_n_flk = min(t_snow_n_flk, tpl_t_f)
773 t_ice_n_flk = min(t_ice_n_flk, tpl_t_f)
774 
775 !_tmp
776 ! Security, avoid too low values (these constraints are used for debugging purposes)
777 ! not good must be deleteted !!!!!!! But the snow scheme is not robust !!!!!!!
778 ! not good must be deleteted !!!!!!! But the snow scheme is not robust !!!!!!!
779 ! not good must be deleteted !!!!!!! But the snow scheme is not robust !!!!!!!
780  t_snow_n_flk = max(t_snow_n_flk, 200.15)
781  t_ice_n_flk = max(t_ice_n_flk, 200.15)
782 !_tmp
783 
784 ! Remove too thin ice and/or snow
785 IF(h_ice_n_flk.LT.h_ice_min_flk) THEN ! Check ice
786  h_ice_n_flk = 0. ! Ice is too thin, remove it, and
787  t_ice_n_flk = tpl_t_f ! set T_ice to the freezing point.
788  h_snow_n_flk = 0. ! Remove snow when there is no ice, and
789  t_snow_n_flk = tpl_t_f ! set T_snow to the freezing point.
790  l_ice_create = .false. ! "Exotic" case, ice has been created but proved to be too thin
791 ELSE IF(h_snow_n_flk.LT.h_snow_min_flk) THEN ! Ice exists, check snow
792  h_snow_n_flk = 0. ! Snow is too thin, remove it,
793  t_snow_n_flk = t_ice_n_flk ! and set the snow temperature equal to the ice temperature.
794 END IF
795 
796 
797 !------------------------------------------------------------------------------
798 ! Compute the mean temperature of the water column.
799 !------------------------------------------------------------------------------
800 
801 IF(l_ice_create) q_w_flk = 0. ! Ice has just been created, set Q_w to zero
802 d_t_mnw_dt = (q_w_flk - q_bot_flk + i_w_flk - i_bot_flk)/tpl_rho_w_r/tpl_c_w/depth_w
803 t_mnw_n_flk = t_mnw_p_flk + d_t_mnw_dt*del_time ! Advance T_mnw
804 t_mnw_n_flk = max(t_mnw_n_flk, tpl_t_f) ! Limit T_mnw by the freezing point
805 
806 
807 !------------------------------------------------------------------------------
808 ! Compute the mixed-layer depth, the mixed-layer temperature,
809 ! the bottom temperature and the shape factor
810 ! with respect to the temperature profile in the thermocline.
811 ! Different formulations are used, depending on the regime of mixing.
812 !------------------------------------------------------------------------------
813 
814 htc_water: IF(h_ice_n_flk.GE.h_ice_min_flk) THEN ! Ice exists
815 
816  t_mnw_n_flk = min(t_mnw_n_flk, tpl_t_r) ! Limit the mean temperature under the ice by T_r
817  t_wml_n_flk = tpl_t_f ! The mixed-layer temperature is equal to the freezing point
818 
819  IF(l_ice_create) THEN ! Ice has just been created
820  IF(h_ml_p_flk.GE.depth_w-h_ml_min_flk) THEN ! h_ML=D when ice is created
821  h_ml_n_flk = 0. ! Set h_ML to zero
822  c_t_n_flk = c_t_min ! Set C_T to its minimum value
823  ELSE ! h_ML<D when ice is created
824  h_ml_n_flk = h_ml_p_flk ! h_ML remains unchanged
825  c_t_n_flk = c_t_p_flk ! C_T (thermocline) remains unchanged
826  END IF
827  t_bot_n_flk = t_wml_n_flk - (t_wml_n_flk-t_mnw_n_flk)/c_t_n_flk/(1.-h_ml_n_flk/depth_w)
828  ! Update the bottom temperature
829 
830  ELSE IF(t_bot_p_flk.LT.tpl_t_r) THEN ! Ice exists and T_bot < T_r, molecular heat transfer
831  h_ml_n_flk = h_ml_p_flk ! h_ML remains unchanged
832  c_t_n_flk = c_t_p_flk ! C_T (thermocline) remains unchanged
833  t_bot_n_flk = t_wml_n_flk - (t_wml_n_flk-t_mnw_n_flk)/c_t_n_flk/(1.-h_ml_n_flk/depth_w)
834  ! Update the bottom temperature
835 
836  ELSE ! Ice exists and T_bot = T_r, convection due to bottom heating
837  t_bot_n_flk = tpl_t_r ! T_bot is equal to the temperature of maximum density
838  IF(h_ml_p_flk.GE.c_small_flk) THEN ! h_ML > 0
839  c_t_n_flk = c_t_p_flk ! C_T (thermocline) remains unchanged
840  h_ml_n_flk = depth_w*(1.-(t_wml_n_flk-t_mnw_n_flk)/(t_wml_n_flk-t_bot_n_flk)/c_t_n_flk)
841  h_ml_n_flk = max(h_ml_n_flk, 0.) ! Update the mixed-layer depth
842  ELSE ! h_ML = 0
843  h_ml_n_flk = h_ml_p_flk ! h_ML remains unchanged
844  c_t_n_flk = (t_wml_n_flk-t_mnw_n_flk)/(t_wml_n_flk-t_bot_n_flk)
845  c_t_n_flk = min(c_t_max, max(c_t_n_flk, c_t_min)) ! Update the shape factor (thermocline)
846  END IF
847  END IF
848 
849  t_bot_n_flk = min(t_bot_n_flk, tpl_t_r) ! Security, limit the bottom temperature by T_r
850 
851 ELSE htc_water ! Open water
852 
853 ! Generalised buoyancy flux scale and convective velocity scale
854  flk_str_1 = flake_buoypar(t_wml_p_flk)*q_star_flk/tpl_rho_w_r/tpl_c_w
855  IF(flk_str_1.LT.0.) THEN
856  w_star_sfc_flk = (-flk_str_1*h_ml_p_flk)**(1./3.) ! Convection
857  ELSE
858  w_star_sfc_flk = 0. ! Neutral or stable stratification
859  END IF
860 
861 !_dm
862 ! The equilibrium depth of the CBL due to surface cooling with the volumetric heating
863 ! is not computed as a solution to the transcendental equation.
864 ! Instead, an algebraic formula is used
865 ! that interpolates between the two asymptotic limits.
866 !_dm
867  conv_equil_h_scale = -q_w_flk/max(i_w_flk, c_small_flk)
868  IF(conv_equil_h_scale.GT.0. .AND. conv_equil_h_scale.LT.1. &
869  .AND. t_wml_p_flk.GT.tpl_t_r) THEN ! The equilibrium CBL depth scale is only used above T_r
870  conv_equil_h_scale = sqrt(6.*conv_equil_h_scale) &
871  + 2.*conv_equil_h_scale/(1.-conv_equil_h_scale)
872  conv_equil_h_scale = min(depth_w, conv_equil_h_scale/extincoef_water_typ)
873  ELSE
874  conv_equil_h_scale = 0. ! Set the equilibrium CBL depth to zero
875  END IF
876 
877 ! Mean buoyancy frequency in the thermocline
878  n_t_mean = flake_buoypar(0.5*(t_wml_p_flk+t_bot_p_flk))*max(0.,(t_wml_p_flk-t_bot_p_flk))
879  IF(h_ml_p_flk.LE.depth_w-h_ml_min_flk) THEN
880  n_t_mean = sqrt(n_t_mean/(depth_w-h_ml_p_flk)) ! Compute N
881  ELSE
882  n_t_mean = 0. ! h_ML=D, set N to zero
883  END IF
884 
885 ! The rate of change of C_T
886  d_c_t_dt = max(w_star_sfc_flk, u_star_w_flk, u_star_min_flk)**2
887  d_c_t_dt = n_t_mean*(depth_w-h_ml_p_flk)**2 &
888  / c_relax_c/d_c_t_dt ! Relaxation time scale for C_T
889  d_c_t_dt = (c_t_max-c_t_min)/max(d_c_t_dt, c_small_flk) ! Rate-of-change of C_T
890  !rsal:
891  d_c_t_dt = min(d_c_t_dt,0.01/del_time)
892 
893 ! Compute the shape factor and the mixed-layer depth,
894 ! using different formulations for convection and wind mixing
895 
896  c_tt_flk = c_tt_1*c_t_p_flk-c_tt_2 ! C_TT, using C_T at the previous time step
897  c_q_flk = 2.*c_tt_flk/c_t_p_flk ! C_Q using C_T at the previous time step
898 
899  mixing_regime: IF(flk_str_1.LT.0.) THEN ! Convective mixing
900 
901  c_t_n_flk = c_t_p_flk + d_c_t_dt*del_time ! Update C_T, assuming dh_ML/dt>0
902  c_t_n_flk = min(c_t_max, max(c_t_n_flk, c_t_min)) ! Limit C_T
903  d_c_t_dt = (c_t_n_flk-c_t_p_flk)/del_time ! Re-compute dC_T/dt
904 
905  IF(h_ml_p_flk.LE.depth_w-h_ml_min_flk) THEN ! Compute dh_ML/dt
906  IF(h_ml_p_flk.LE.h_ml_min_flk) THEN ! Use a reduced entrainment equation (spin-up)
907  d_h_ml_dt = c_cbl_1/c_cbl_2*max(w_star_sfc_flk, c_small_flk)
908 
909 !_dbg
910 ! WRITE*, ' FLake: reduced entrainment eq. D_time*d_h_ML_dt = ', d_h_ML_dt*del_time
911 ! WRITE*, ' w_* = ', w_star_sfc_flk
912 ! WRITE*, ' \beta*Q_* = ', flk_str_1
913 !_dbg
914 
915  ELSE ! Use a complete entrainment equation
916  r_h_icesnow = depth_w/h_ml_p_flk
917  r_rho_c_icesnow = r_h_icesnow-1.
918  r_ti_icesnow = c_t_p_flk/c_tt_flk
919  r_tstar_icesnow = (r_ti_icesnow/2.-1.)*r_rho_c_icesnow + 1.
920  d_h_ml_dt = -q_star_flk*(r_tstar_icesnow*(1.+c_cbl_1)-1.) - q_bot_flk
921  d_h_ml_dt = d_h_ml_dt/tpl_rho_w_r/tpl_c_w ! Q_* and Q_b flux terms
922  flk_str_2 = (depth_w-h_ml_p_flk)*(t_wml_p_flk-t_bot_p_flk)*c_tt_2/c_tt_flk*d_c_t_dt
923  d_h_ml_dt = d_h_ml_dt + flk_str_2 ! Add dC_T/dt term
924  flk_str_2 = i_bot_flk + (r_ti_icesnow-1.)*i_h_flk - r_ti_icesnow*i_intm_h_d_flk
925  flk_str_2 = flk_str_2 + (r_ti_icesnow-2.)*r_rho_c_icesnow*(i_h_flk-i_intm_0_h_flk)
926  flk_str_2 = flk_str_2/tpl_rho_w_r/tpl_c_w
927  d_h_ml_dt = d_h_ml_dt + flk_str_2 ! Add radiation terms
928  flk_str_2 = -c_cbl_2*r_tstar_icesnow*q_star_flk/tpl_rho_w_r/tpl_c_w/max(w_star_sfc_flk, c_small_flk)
929  flk_str_2 = flk_str_2 + c_t_p_flk*(t_wml_p_flk-t_bot_p_flk)
930  d_h_ml_dt = d_h_ml_dt/flk_str_2 ! dh_ML/dt = r.h.s.
931  END IF
932 !_dm
933 ! Notice that dh_ML/dt may appear to be negative
934 ! (e.g. due to buoyancy loss to bottom sediments and/or
935 ! the effect of volumetric radiation heating),
936 ! although a negative generalized buoyancy flux scale indicates
937 ! that the equilibrium CBL depth has not yet been reached
938 ! and convective deepening of the mixed layer should take place.
939 ! Physically, this situation reflects an approximate character of the lake model.
940 ! Using the self-similar temperature profile in the thermocline,
941 ! there is always communication between the mixed layer, the thermocline
942 ! and the lake bottom. As a result, the rate of change of the CBL depth
943 ! is always dependent on the bottom heat flux and the radiation heating of the thermocline.
944 ! In reality, convective mixed-layer deepening may be completely decoupled
945 ! from the processes underneath. In order to account for this fact,
946 ! the rate of CBL deepening is set to a small value
947 ! if dh_ML/dt proves to be negative.
948 ! This is "double insurance" however,
949 ! as a negative dh_ML/dt is encountered very rarely.
950 !_dm
951 
952  d_h_ml_dt = max(d_h_ml_dt, c_small_flk)
953  h_ml_n_flk = h_ml_p_flk + d_h_ml_dt*del_time ! Update h_ML
954  h_ml_n_flk = max(h_ml_min_flk, min(h_ml_n_flk, depth_w)) ! Security, limit h_ML
955  ELSE ! Mixing down to the lake bottom
956  h_ml_n_flk = depth_w
957  END IF
958 
959  ELSE mixing_regime ! Wind mixing
960 
961  d_h_ml_dt = max(u_star_w_flk, u_star_min_flk) ! The surface friction velocity
962  zm_h_scale = (abs(par_coriolis)/c_sbl_zm_n + n_t_mean/c_sbl_zm_i)*d_h_ml_dt**2
963  zm_h_scale = zm_h_scale + flk_str_1/c_sbl_zm_s
964  zm_h_scale = max(zm_h_scale, c_small_flk)
965  zm_h_scale = d_h_ml_dt**3/zm_h_scale
966  zm_h_scale = max(h_ml_min_flk, min(zm_h_scale, h_ml_max_flk)) ! The ZM96 SBL depth scale
967  zm_h_scale = max(zm_h_scale, conv_equil_h_scale) ! Equilibrium mixed-layer depth
968 
969 !_dm
970 ! In order to avoid numerical discretization problems,
971 ! an analytical solution to the evolution equation
972 ! for the wind-mixed layer depth is used.
973 ! That is, an exponential relaxation formula is applied
974 ! over the time interval equal to the model time step.
975 !_dm
976 
977  d_h_ml_dt = c_relax_h*d_h_ml_dt/zm_h_scale*del_time
978  h_ml_n_flk = zm_h_scale - (zm_h_scale-h_ml_p_flk)*exp(-d_h_ml_dt) ! Update h_ML
979  h_ml_n_flk = max(h_ml_min_flk, min(h_ml_n_flk, depth_w)) ! Limit h_ML
980  d_h_ml_dt = (h_ml_n_flk-h_ml_p_flk)/del_time ! Re-compute dh_ML/dt
981 
982  IF(h_ml_n_flk.LE.h_ml_p_flk) &
983  d_c_t_dt = -d_c_t_dt ! Mixed-layer retreat or stationary state, dC_T/dt<0
984  c_t_n_flk = c_t_p_flk + d_c_t_dt*del_time ! Update C_T
985  c_t_n_flk = min(c_t_max, max(c_t_n_flk, c_t_min)) ! Limit C_T
986  d_c_t_dt = (c_t_n_flk-c_t_p_flk)/del_time ! Re-compute dC_T/dt
987 
988 
989  END IF mixing_regime
990 
991 ! Compute the time-rate-of-change of the the bottom temperature,
992 ! depending on the sign of dh_ML/dt
993 ! Update the bottom temperature and the mixed-layer temperature
994 
995  IF(h_ml_n_flk.LE.depth_w-h_ml_min_flk) THEN ! Mixing did not reach the bottom
996 
997  IF(h_ml_n_flk.GT.h_ml_p_flk) THEN ! Mixed-layer deepening
998  r_h_icesnow = h_ml_p_flk/depth_w
999  r_rho_c_icesnow = 1.-r_h_icesnow
1000  r_ti_icesnow = 0.5*c_t_p_flk*r_rho_c_icesnow+c_tt_flk*(2.*r_h_icesnow-1.)
1001  r_tstar_icesnow = (0.5+c_tt_flk-c_q_flk)/r_ti_icesnow
1002  r_ti_icesnow = (1.-c_t_p_flk*r_rho_c_icesnow)/r_ti_icesnow
1003 
1004  d_t_bot_dt = (q_w_flk-q_bot_flk+i_w_flk-i_bot_flk)/tpl_rho_w_r/tpl_c_w
1005  d_t_bot_dt = d_t_bot_dt - c_t_p_flk*(t_wml_p_flk-t_bot_p_flk)*d_h_ml_dt
1006  d_t_bot_dt = d_t_bot_dt*r_tstar_icesnow/depth_w ! Q+I fluxes and dh_ML/dt term
1007 
1008  flk_str_2 = i_intm_h_d_flk - (1.-c_q_flk)*i_h_flk - c_q_flk*i_bot_flk
1009  flk_str_2 = flk_str_2*r_ti_icesnow/(depth_w-h_ml_p_flk)/tpl_rho_w_r/tpl_c_w
1010  d_t_bot_dt = d_t_bot_dt + flk_str_2 ! Add radiation-flux term
1011 
1012  flk_str_2 = (1.-c_tt_2*r_ti_icesnow)/c_t_p_flk
1013  flk_str_2 = flk_str_2*(t_wml_p_flk-t_bot_p_flk)*d_c_t_dt
1014  d_t_bot_dt = d_t_bot_dt + flk_str_2 ! Add dC_T/dt term
1015 
1016  ELSE ! Mixed-layer retreat or stationary state
1017  d_t_bot_dt = 0. ! dT_bot/dt=0
1018  END IF
1019 
1020  t_bot_n_flk = t_bot_p_flk + d_t_bot_dt*del_time ! Update T_bot
1021  t_bot_n_flk = max(t_bot_n_flk, tpl_t_f) ! Security, limit T_bot by the freezing point
1022  flk_str_2 = (t_bot_n_flk-tpl_t_r)*flake_buoypar(t_mnw_n_flk)
1023  IF(flk_str_2.LT.0.) t_bot_n_flk = tpl_t_r ! Security, avoid T_r crossover
1024  t_wml_n_flk = c_t_n_flk*(1.-h_ml_n_flk/depth_w)
1025  t_wml_n_flk = (t_mnw_n_flk-t_bot_n_flk*t_wml_n_flk)/(1.-t_wml_n_flk)
1026  t_wml_n_flk = max(t_wml_n_flk, tpl_t_f) ! Security, limit T_wML by the freezing point
1027 
1028  ELSE ! Mixing down to the lake bottom
1029 
1030  h_ml_n_flk = depth_w
1031  t_wml_n_flk = t_mnw_n_flk
1032  t_bot_n_flk = t_mnw_n_flk
1033  c_t_n_flk = c_t_min
1034 !rsal:
1035  c_t_n_flk = 0.65
1036 
1037 
1038  END IF
1039 
1040 END IF htc_water
1041 
1042 !------------------------------------------------------------------------------
1043 ! Compute the depth of the upper layer of bottom sediments
1044 ! and the temperature at that depth.
1045 !------------------------------------------------------------------------------
1046 
1047 !SURFEX
1048 sediment: IF(lflk_botsed_use) THEN ! The bottom-sediment scheme is used
1049 
1050  IF(h_b1_p_flk.GE.depth_bs-h_b1_min_flk) THEN ! No T(z) maximum (no thermal wave)
1051  h_b1_p_flk = 0. ! Set H_B1_p to zero
1052  t_b1_p_flk = t_bot_p_flk ! Set T_B1_p to the bottom temperature
1053  END IF
1054 
1055  flk_str_1 = 2.*phi_b1_pr0/(1.-c_b1)*tpl_kappa_w/tpl_rho_w_r/tpl_c_w*del_time
1056  h_ice_threshold = sqrt(flk_str_1) ! Threshold value of H_B1
1057  h_ice_threshold = min(0.9*depth_bs, h_ice_threshold) ! Limit H_B1
1058  flk_str_2 = c_b2/(1.-c_b2)*(t_bs-t_b1_p_flk)/(depth_bs-h_b1_p_flk)
1059 
1060  IF(h_b1_p_flk.LT.h_ice_threshold) THEN ! Use a truncated equation for H_B1(t)
1061  h_b1_n_flk = sqrt(h_b1_p_flk**2+flk_str_1) ! Advance H_B1
1062  d_h_b1_dt = (h_b1_n_flk-h_b1_p_flk)/del_time ! Re-compute dH_B1/dt
1063  ELSE ! Use a full equation for H_B1(t)
1064  flk_str_1 = (q_bot_flk+i_bot_flk)/h_b1_p_flk/tpl_rho_w_r/tpl_c_w
1065  flk_str_1 = flk_str_1 - (1.-c_b1)*(t_bot_n_flk-t_bot_p_flk)/del_time
1066  d_h_b1_dt = (1.-c_b1)*(t_bot_p_flk-t_b1_p_flk)/h_b1_p_flk + c_b1*flk_str_2
1067  d_h_b1_dt = flk_str_1/d_h_b1_dt
1068  h_b1_n_flk = h_b1_p_flk + d_h_b1_dt*del_time ! Advance H_B1
1069  END IF
1070  d_t_b1_dt = flk_str_2*d_h_b1_dt
1071  t_b1_n_flk = t_b1_p_flk + d_t_b1_dt*del_time ! Advance T_B1
1072 
1073 ! Use a very simplistic procedure, where only the upper layer profile is used,
1074 ! H_B1 is always set to depth_bs, and T_B1 is always set to T_bs.
1075 ! Then, the time derivatives are zero, and the sign of the bottom heat flux depends on
1076 ! whether T_bot is smaller or greater than T_bs.
1077 ! This is, of course, an oversimplified scheme.
1078 
1079  l_snow_exists = h_b1_n_flk.GE.depth_bs-h_b1_min_flk &! H_B1 reached depth_bs, or
1080  .OR. h_b1_n_flk.LT.h_b1_min_flk &! H_B1 decreased to zero, or
1081  .OR.(t_bot_n_flk-t_b1_n_flk)*(t_bs-t_b1_n_flk).LE.0. ! there is no T(z) maximum
1082  IF(l_snow_exists) THEN
1083  h_b1_n_flk = depth_bs ! Set H_B1 to the depth of the thermally active layer
1084  t_b1_n_flk = t_bs ! Set T_B1 to the climatological temperature
1085  END IF
1086 
1087 ELSE sediment ! The bottom-sediment scheme is not used
1088 
1089  h_b1_n_flk = rflk_depth_bs_ref ! H_B1 is set to a reference value
1090  t_b1_n_flk = tpl_t_r ! T_B1 is set to the temperature of maximum density
1091 
1092 END IF sediment
1093 
1094 
1095 !------------------------------------------------------------------------------
1096 ! Impose additional constraints.
1097 !------------------------------------------------------------------------------
1098 
1099 ! In case of unstable stratification, force mixing down to the bottom
1100 flk_str_2 = (t_wml_n_flk-t_bot_n_flk)*flake_buoypar(t_mnw_n_flk)
1101 IF(flk_str_2.LT.0.) THEN
1102 
1103  h_ml_n_flk = depth_w
1104  t_wml_n_flk = t_mnw_n_flk
1105  t_bot_n_flk = t_mnw_n_flk
1106  c_t_n_flk = c_t_min
1107 !rsal:
1108  c_t_n_flk = 0.65
1109 
1110 END IF
1111 
1112 
1113 !------------------------------------------------------------------------------
1114 ! Update the surface temperature.
1115 !------------------------------------------------------------------------------
1116 
1117 IF(h_snow_n_flk.GE.h_snow_min_flk) THEN
1118  t_sfc_n = t_snow_n_flk ! Snow exists, use the snow temperature
1119 ELSE IF(h_ice_n_flk.GE.h_ice_min_flk) THEN
1120  t_sfc_n = t_ice_n_flk ! Ice exists but there is no snow, use the ice temperature
1121 ELSE
1122  t_sfc_n = t_wml_n_flk ! No ice-snow cover, use the mixed-layer temperature
1123 END IF
1124 IF (lhook) CALL dr_hook('FLAKE:FLAKE_DRIVER',1,zhook_handle)
1125 
1126 !------------------------------------------------------------------------------
1127 ! End calculations
1128 !==============================================================================
1129 
1130 END SUBROUTINE flake_driver
1131 
1132 
1133 !==============================================================================
1134 
1135 !==============================================================================
1136 !SURFEX include 'flake_buoypar.incf'
1137 ! File %M% from Library %Q%
1138 ! Version %I% from %G% extracted: %H%
1139 !------------------------------------------------------------------------------
1140 
1141 !SURFEX REAL FUNCTION flake_buoypar (T_water)
1142 FUNCTION flake_buoypar (T_water)
1143 
1144 !------------------------------------------------------------------------------
1145 !
1146 ! Description:
1147 !
1148 ! Computes the buoyancy parameter,
1149 ! using a quadratic equation of state for the fresh-water.
1150 !
1151 !
1152 ! Current Code Owner: DWD, Dmitrii Mironov
1153 ! Phone: +49-69-8062 2705
1154 ! Fax: +49-69-8062 3721
1155 ! E-mail: dmitrii.mironov@dwd.de
1156 !
1157 ! History:
1158 ! Version Date Name
1159 ! ---------- ---------- ----
1160 ! 1.00 2005/11/17 Dmitrii Mironov
1161 ! Initial release
1162 ! !VERSION! !DATE! <Your name>
1163 ! <Modification comments>
1164 !
1165 ! Code Description:
1166 ! Language: Fortran 90.
1167 ! Software Standards: "European Standards for Writing and
1168 ! Documenting Exchangeable Fortran 90 Code".
1169 !==============================================================================
1170 !
1171 ! Declarations:
1172 !
1173 ! Modules used:
1174 
1175 !_dm Parameters are USEd in module "flake".
1176 !_nu USE modd_data_parameters , ONLY : &
1177 !_nu ireals, & ! KIND-type parameter for real variables
1178 !_nu iintegers ! KIND-type parameter for "normal" integer variables
1179 
1180 USE modd_flake_parameters , ONLY : &
1181  tpl_grav , &! Acceleration due to gravity [m s^{-2}]
1182  tpl_t_r , &! Temperature of maximum density of fresh water [K]
1183  tpl_a_t ! Constant in the fresh-water equation of state [K^{-2}]
1184 
1185 !==============================================================================
1186 
1187 IMPLICIT NONE
1188 
1189 !==============================================================================
1190 !
1191 ! Declarations
1192 
1193 ! Input (function argument)
1194 REAL , INTENT(IN) :: &
1195  t_water ! Water temperature [K]
1196 
1197 ! Output (function result)
1198 REAL :: &
1199  flake_buoypar ! Buoyancy parameter [m s^{-2} K^{-1}]
1200 REAL(KIND=JPRB) :: zhook_handle
1201 
1202 !==============================================================================
1203 ! Start calculations
1204 !------------------------------------------------------------------------------
1205 
1206 ! Buoyancy parameter [m s^{-2} K^{-1}]
1207 
1208  IF (lhook) CALL dr_hook('FLAKE:FLAKE_BUOYPAR',0,zhook_handle)
1209  flake_buoypar = tpl_grav*tpl_a_t*(t_water-tpl_t_r)
1210 IF (lhook) CALL dr_hook('FLAKE:FLAKE_BUOYPAR',1,zhook_handle)
1211 
1212 !------------------------------------------------------------------------------
1213 ! End calculations
1214 !==============================================================================
1215 
1216 END FUNCTION flake_buoypar
1217 
1218 
1219 !==============================================================================
1220 !SURFEX include 'flake_snowdensity.incf'
1221 ! File %M% from Library %Q%
1222 ! Version %I% from %G% extracted: %H%
1223 !------------------------------------------------------------------------------
1224 
1225 !SURFEX REAL FUNCTION flake_snowdensity (h_snow)
1226 FUNCTION flake_snowdensity (h_snow)
1227 
1228 !------------------------------------------------------------------------------
1229 !
1230 ! Description:
1231 !
1232 ! Computes the snow density,
1233 ! using an empirical approximation from Heise et al. (2003).
1234 !
1235 !
1236 ! Current Code Owner: DWD, Dmitrii Mironov
1237 ! Phone: +49-69-8062 2705
1238 ! Fax: +49-69-8062 3721
1239 ! E-mail: dmitrii.mironov@dwd.de
1240 !
1241 ! History:
1242 ! Version Date Name
1243 ! ---------- ---------- ----
1244 ! 1.00 2005/11/17 Dmitrii Mironov
1245 ! Initial release
1246 ! !VERSION! !DATE! <Your name>
1247 ! <Modification comments>
1248 !
1249 ! Code Description:
1250 ! Language: Fortran 90.
1251 ! Software Standards: "European Standards for Writing and
1252 ! Documenting Exchangeable Fortran 90 Code".
1253 !==============================================================================
1254 !
1255 ! Declarations:
1256 !
1257 ! Modules used:
1258 
1259 !_dm Parameters are USEd in module "flake".
1260 !_nu USE modd_data_parameters , ONLY : &
1261 !_nu ireals, & ! KIND-type parameter for real variables
1262 !_nu iintegers ! KIND-type parameter for "normal" integer variables
1263 
1264 USE modd_flake_parameters , ONLY : &
1265  tpl_rho_w_r , &! Maximum density of fresh water [kg m^{-3}]
1266  tpl_rho_s_min , &! Minimum snow density [kg m^{-3}]
1267  tpl_rho_s_max , &! Maximum snow density [kg m^{-3}]
1268  tpl_gamma_rho_s , &! Empirical parameter [kg m^{-4}] in the expression for the snow density
1269  c_small_flk ! A small number
1270 
1271 !==============================================================================
1272 
1273 IMPLICIT NONE
1274 
1275 !==============================================================================
1276 !
1277 ! Declarations
1278 
1279 ! Input (function argument)
1280 REAL , INTENT(IN) :: &
1281  h_snow ! Snow thickness [m]
1282 
1283 ! Output (function result)
1284 REAL :: &
1285  flake_snowdensity ! Snow density [kg m^{-3}]
1286 REAL(KIND=JPRB) :: zhook_handle
1287 
1288 !==============================================================================
1289 ! Start calculations
1290 !------------------------------------------------------------------------------
1291 
1292 ! Snow density [kg m^{-3}]
1293 
1294 ! Security. Ensure that the expression in () does not become negative at a very large h_snow.
1295  IF (lhook) CALL dr_hook('FLAKE:FLAKE_SNOWDENSITY',0,zhook_handle)
1296  flake_snowdensity = max( c_small_flk, (1. - h_snow*tpl_gamma_rho_s/tpl_rho_w_r) )
1297  flake_snowdensity = min( tpl_rho_s_max, tpl_rho_s_min/flake_snowdensity )
1298 IF (lhook) CALL dr_hook('FLAKE:FLAKE_SNOWDENSITY',1,zhook_handle)
1299 
1300 !------------------------------------------------------------------------------
1301 ! End calculations
1302 !==============================================================================
1303 
1304 END FUNCTION flake_snowdensity
1305 
1306 
1307 !==============================================================================
1308 !SURFEX include 'flake_snowheatconduct.incf'
1309 ! File %M% from Library %Q%
1310 ! Version %I% from %G% extracted: %H%
1311 !------------------------------------------------------------------------------
1312 
1313 !SURFEX REAL FUNCTION flake_snowheatconduct (h_snow)
1314 FUNCTION flake_snowheatconduct (h_snow)
1315 
1316 !------------------------------------------------------------------------------
1317 !
1318 ! Description:
1319 !
1320 ! Computes the snow heat conductivity,
1321 ! using an empirical approximation from Heise et al. (2003).
1322 !
1323 !
1324 ! Current Code Owner: DWD, Dmitrii Mironov
1325 ! Phone: +49-69-8062 2705
1326 ! Fax: +49-69-8062 3721
1327 ! E-mail: dmitrii.mironov@dwd.de
1328 !
1329 ! History:
1330 ! Version Date Name
1331 ! ---------- ---------- ----
1332 ! 1.00 2005/11/17 Dmitrii Mironov
1333 ! Initial release
1334 ! !VERSION! !DATE! <Your name>
1335 ! <Modification comments>
1336 !
1337 ! Code Description:
1338 ! Language: Fortran 90.
1339 ! Software Standards: "European Standards for Writing and
1340 ! Documenting Exchangeable Fortran 90 Code".
1341 !==============================================================================
1342 !
1343 ! Declarations:
1344 !
1345 ! Modules used:
1346 
1347 !_dm Parameters are USEd in module "flake".
1348 !_nu USE modd_data_parameters , ONLY : &
1349 !_nu ireals, & ! KIND-type parameter for real variables
1350 !_nu iintegers ! KIND-type parameter for "normal" integer variables
1351 
1352 USE modd_flake_parameters , ONLY : &
1353  tpl_rho_w_r , &! Maximum density of fresh water [kg m^{-3}]
1354  tpl_kappa_s_min , &! Minimum molecular heat conductivity of snow [J m^{-1} s^{-1} K^{-1}]
1355  tpl_kappa_s_max , &! Maximum molecular heat conductivity of snow [J m^{-1} s^{-1} K^{-1}]
1356  tpl_gamma_kappa_s ! Empirical parameter [J m^{-2} s^{-1} K^{-1}] in the expression for kappa_S
1357 
1358 !==============================================================================
1359 
1360 IMPLICIT NONE
1361 
1362 !==============================================================================
1363 !
1364 ! Declarations
1365 
1366 ! Input (function argument)
1367 REAL , INTENT(IN) :: &
1368  h_snow ! Snow thickness [m]
1369 
1370 ! Output (function result)
1371 REAL :: &
1372  flake_snowheatconduct ! Snow heat conductivity [J m^{-1} s^{-1} K^{-1} = kg m s^{-3} K^{-1}]
1373 REAL(KIND=JPRB) :: zhook_handle
1374 
1375 !==============================================================================
1376 ! Start calculations
1377 !------------------------------------------------------------------------------
1378 
1379 ! Snow heat conductivity [J m^{-1} s^{-1} K^{-1} = kg m s^{-3} K^{-1}]
1380 
1381  IF (lhook) CALL dr_hook('FLAKE:FLAKE_SNOWHEATCONDUCT',0,zhook_handle)
1382  flake_snowheatconduct = flake_snowdensity( h_snow ) ! Compute snow density
1383  flake_snowheatconduct = min( tpl_kappa_s_max, tpl_kappa_s_min &
1384  + h_snow*tpl_gamma_kappa_s*flake_snowheatconduct/tpl_rho_w_r )
1385 IF (lhook) CALL dr_hook('FLAKE:FLAKE_SNOWHEATCONDUCT',1,zhook_handle)
1386 
1387 !------------------------------------------------------------------------------
1388 ! End calculations
1389 !==============================================================================
1390 
1391 END FUNCTION flake_snowheatconduct
1392 
1393 END MODULE mode_flake
1394 
subroutine flake_radflux(depth_w, albedo, opticpar_water, opticpar_ice, opticpar_snow)
Definition: mode_flake.F90:237
real function flake_snowdensity(h_snow)
real function flake_buoypar(T_water)
real function flake_snowheatconduct(h_snow)
subroutine flake_driver(depth_w, depth_bs, T_bs, par_Coriolis, extincoef_water_typ, del_time, T_sfc_p, T_sfc_n)
Definition: mode_flake.F90:399