source: palm/trunk/SOURCE/microphysics_mod.f90 @ 2000

Last change on this file since 2000 was 2000, checked in by knoop, 8 years ago

Forced header and separation lines into 80 columns

  • Property svn:keywords set to Id
File size: 86.7 KB
Line 
1!> @file microphysics_mod.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
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2016 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22! Forced header and separation lines into 80 columns
23!
24! Former revisions:
25! -----------------
26! $Id: microphysics_mod.f90 2000 2016-08-20 18:09:15Z knoop $
27!
28! 1850 2016-04-08 13:29:27Z maronga
29! Module renamed
30! Adapted for modularization of microphysics.
31!
32! 1845 2016-04-08 08:29:13Z raasch
33! nzb_2d replaced by nzb_s_inner, Kessler precipitation is stored at surface
34! point (instead of one point above surface)
35!
36! 1831 2016-04-07 13:15:51Z hoffmann
37! turbulence renamed collision_turbulence,
38! drizzle renamed cloud_water_sedimentation. cloud_water_sedimentation also
39! avaialble for microphysics_kessler.
40!
41! 1822 2016-04-07 07:49:42Z hoffmann
42! Unused variables removed.
43! Kessler scheme integrated.
44!
45! 1691 2015-10-26 16:17:44Z maronga
46! Added new routine calc_precipitation_amount. The routine now allows to account
47! for precipitation due to sedimenation of cloud (fog) droplets
48!
49! 1682 2015-10-07 23:56:08Z knoop
50! Code annotations made doxygen readable
51!
52! 1646 2015-09-02 16:00:10Z hoffmann
53! Bugfix: Wrong computation of d_mean.
54!
55! 1361 2014-04-16 15:17:48Z hoffmann
56! Bugfix in sedimentation_rain: Index corrected.
57! Vectorized version of adjust_cloud added.
58! Little reformatting of the code.
59!
60! 1353 2014-04-08 15:21:23Z heinze
61! REAL constants provided with KIND-attribute
62!
63! 1346 2014-03-27 13:18:20Z heinze
64! Bugfix: REAL constants provided with KIND-attribute especially in call of
65! intrinsic function like MAX, MIN, SIGN
66!
67! 1334 2014-03-25 12:21:40Z heinze
68! Bugfix: REAL constants provided with KIND-attribute
69!
70! 1322 2014-03-20 16:38:49Z raasch
71! REAL constants defined as wp-kind
72!
73! 1320 2014-03-20 08:40:49Z raasch
74! ONLY-attribute added to USE-statements,
75! kind-parameters added to all INTEGER and REAL declaration statements,
76! kinds are defined in new module kinds,
77! comment fields (!:) to be used for variable explanations added to
78! all variable declaration statements
79!
80! 1241 2013-10-30 11:36:58Z heinze
81! hyp and rho have to be calculated at each time step if data from external
82! file LSF_DATA are used
83!
84! 1115 2013-03-26 18:16:16Z hoffmann
85! microphyical tendencies are calculated in microphysics_control in an optimized
86! way; unrealistic values are prevented; bugfix in evaporation; some reformatting
87!
88! 1106 2013-03-04 05:31:38Z raasch
89! small changes in code formatting
90!
91! 1092 2013-02-02 11:24:22Z raasch
92! unused variables removed
93! file put under GPL
94!
95! 1065 2012-11-22 17:42:36Z hoffmann
96! Sedimentation process implemented according to Stevens and Seifert (2008).
97! Turbulence effects on autoconversion and accretion added (Seifert, Nuijens
98! and Stevens, 2010).
99!
100! 1053 2012-11-13 17:11:03Z hoffmann
101! initial revision
102!
103! Description:
104! ------------
105!> Calculate bilk cloud microphysics.
106!------------------------------------------------------------------------------!
107 MODULE microphysics_mod
108
109    USE  kinds
110
111    IMPLICIT NONE
112
113    LOGICAL ::  cloud_water_sedimentation = .FALSE.  !< cloud water sedimentation
114    LOGICAL ::  limiter_sedimentation = .TRUE.       !< sedimentation limiter
115    LOGICAL ::  collision_turbulence = .FALSE.       !< turbulence effects
116    LOGICAL ::  ventilation_effect = .TRUE.          !< ventilation effect
117
118    REAL(wp) ::  a_1 = 8.69E-4_wp          !< coef. in turb. parametrization (cm-2 s3)
119    REAL(wp) ::  a_2 = -7.38E-5_wp         !< coef. in turb. parametrization (cm-2 s3)
120    REAL(wp) ::  a_3 = -1.40E-2_wp         !< coef. in turb. parametrization
121    REAL(wp) ::  a_term = 9.65_wp          !< coef. for terminal velocity (m s-1)
122    REAL(wp) ::  a_vent = 0.78_wp          !< coef. for ventilation effect
123    REAL(wp) ::  b_1 = 11.45E-6_wp         !< coef. in turb. parametrization (m)
124    REAL(wp) ::  b_2 = 9.68E-6_wp          !< coef. in turb. parametrization (m)
125    REAL(wp) ::  b_3 = 0.62_wp             !< coef. in turb. parametrization
126    REAL(wp) ::  b_term = 9.8_wp           !< coef. for terminal velocity (m s-1)
127    REAL(wp) ::  b_vent = 0.308_wp         !< coef. for ventilation effect
128    REAL(wp) ::  beta_cc = 3.09E-4_wp      !< coef. in turb. parametrization (cm-2 s3)
129    REAL(wp) ::  c_1 = 4.82E-6_wp          !< coef. in turb. parametrization (m)
130    REAL(wp) ::  c_2 = 4.8E-6_wp           !< coef. in turb. parametrization (m)
131    REAL(wp) ::  c_3 = 0.76_wp             !< coef. in turb. parametrization
132    REAL(wp) ::  c_const = 0.93_wp         !< const. in Taylor-microscale Reynolds number
133    REAL(wp) ::  c_evap = 0.7_wp           !< constant in evaporation
134    REAL(wp) ::  c_term = 600.0_wp         !< coef. for terminal velocity (m-1)
135    REAL(wp) ::  diff_coeff_l = 0.23E-4_wp  !< diffusivity of water vapor (m2 s-1)
136    REAL(wp) ::  eps_sb = 1.0E-20_wp       !< threshold in two-moments scheme
137    REAL(wp) ::  k_cc = 9.44E09_wp         !< const. cloud-cloud kernel (m3 kg-2 s-1)
138    REAL(wp) ::  k_cr0 = 4.33_wp           !< const. cloud-rain kernel (m3 kg-1 s-1)
139    REAL(wp) ::  k_rr = 7.12_wp            !< const. rain-rain kernel (m3 kg-1 s-1)
140    REAL(wp) ::  k_br = 1000.0_wp          !< const. in breakup parametrization (m-1)
141    REAL(wp) ::  k_st = 1.2E8_wp           !< const. in drizzle parametrization (m-1 s-1)
142    REAL(wp) ::  kappa_rr = 60.7_wp        !< const. in collision kernel (kg-1/3)
143    REAL(wp) ::  kin_vis_air = 1.4086E-5_wp  !< kin. viscosity of air (m2 s-1)
144    REAL(wp) ::  prec_time_const = 0.001_wp  !< coef. in Kessler scheme (s-1)
145    REAL(wp) ::  ql_crit = 0.0005_wp       !< coef. in Kessler scheme (kg kg-1)
146    REAL(wp) ::  schmidt_p_1d3=0.8921121_wp  !< Schmidt number**0.33333, 0.71**0.33333
147    REAL(wp) ::  sigma_gc = 1.3_wp         !< geometric standard deviation cloud droplets
148    REAL(wp) ::  thermal_conductivity_l = 2.43E-2_wp  !< therm. cond. air (J m-1 s-1 K-1)
149    REAL(wp) ::  w_precipitation = 9.65_wp  !< maximum terminal velocity (m s-1)
150    REAL(wp) ::  x0 = 2.6E-10_wp           !< separating drop mass (kg)
151    REAL(wp) ::  xrmin = 2.6E-10_wp        !< minimum rain drop size (kg)
152    REAL(wp) ::  xrmax = 5.0E-6_wp         !< maximum rain drop site (kg)
153
154    REAL(wp) ::  c_sedimentation = 2.0_wp  !< Courant number of sedimentation process
155    REAL(wp) ::  dpirho_l                  !< 6.0 / ( pi * rho_l )
156    REAL(wp) ::  dt_micro                  !< microphysics time step
157    REAL(wp) ::  nc_const = 70.0E6_wp      !< cloud droplet concentration
158    REAL(wp) ::  dt_precipitation = 100.0_wp !< timestep precipitation (s)
159    REAL(wp) ::  sed_qc_const              !< const. for sedimentation of cloud water
160    REAL(wp) ::  pirho_l                   !< pi * rho_l / 6.0;
161
162    REAL(wp), DIMENSION(:), ALLOCATABLE ::  nc_1d  !<
163    REAL(wp), DIMENSION(:), ALLOCATABLE ::  nr_1d  !<
164    REAL(wp), DIMENSION(:), ALLOCATABLE ::  pt_1d  !<
165    REAL(wp), DIMENSION(:), ALLOCATABLE ::  qc_1d  !<
166    REAL(wp), DIMENSION(:), ALLOCATABLE ::  qr_1d  !<
167    REAL(wp), DIMENSION(:), ALLOCATABLE ::  q_1d   !<
168
169    SAVE
170
171    PRIVATE
172    PUBLIC microphysics_control, microphysics_init
173
174    PUBLIC cloud_water_sedimentation, collision_turbulence, c_sedimentation,   &
175           dt_precipitation, limiter_sedimentation, nc_const, sigma_gc,        &
176           ventilation_effect
177
178    INTERFACE microphysics_control
179       MODULE PROCEDURE microphysics_control
180       MODULE PROCEDURE microphysics_control_ij
181    END INTERFACE microphysics_control
182
183    INTERFACE adjust_cloud
184       MODULE PROCEDURE adjust_cloud
185       MODULE PROCEDURE adjust_cloud_ij
186    END INTERFACE adjust_cloud
187
188    INTERFACE autoconversion
189       MODULE PROCEDURE autoconversion
190       MODULE PROCEDURE autoconversion_ij
191    END INTERFACE autoconversion
192
193    INTERFACE autoconversion_kessler
194       MODULE PROCEDURE autoconversion_kessler
195       MODULE PROCEDURE autoconversion_kessler_ij
196    END INTERFACE autoconversion_kessler
197
198    INTERFACE accretion
199       MODULE PROCEDURE accretion
200       MODULE PROCEDURE accretion_ij
201    END INTERFACE accretion
202
203    INTERFACE selfcollection_breakup
204       MODULE PROCEDURE selfcollection_breakup
205       MODULE PROCEDURE selfcollection_breakup_ij
206    END INTERFACE selfcollection_breakup
207
208    INTERFACE evaporation_rain
209       MODULE PROCEDURE evaporation_rain
210       MODULE PROCEDURE evaporation_rain_ij
211    END INTERFACE evaporation_rain
212
213    INTERFACE sedimentation_cloud
214       MODULE PROCEDURE sedimentation_cloud
215       MODULE PROCEDURE sedimentation_cloud_ij
216    END INTERFACE sedimentation_cloud
217 
218    INTERFACE sedimentation_rain
219       MODULE PROCEDURE sedimentation_rain
220       MODULE PROCEDURE sedimentation_rain_ij
221    END INTERFACE sedimentation_rain
222
223    INTERFACE calc_precipitation_amount
224       MODULE PROCEDURE calc_precipitation_amount
225       MODULE PROCEDURE calc_precipitation_amount_ij
226    END INTERFACE calc_precipitation_amount
227
228 CONTAINS
229!------------------------------------------------------------------------------!
230! Description:
231! ------------
232!> Initialization of bulk microphysics
233!------------------------------------------------------------------------------!
234    SUBROUTINE microphysics_init
235
236       USE arrays_3d,                                                          &
237           ONLY:  dzu
238
239       USE constants,                                                          &
240           ONLY:  pi
241
242       USE cloud_parameters,                                                   &
243           ONLY:  rho_l
244
245       USE control_parameters,                                                 &
246           ONLY:  microphysics_seifert
247
248       USE indices,                                                            &
249           ONLY:  nzb, nzt
250
251       IMPLICIT NONE
252
253!
254!--    constant for the sedimentation of cloud water (2-moment cloud physics)
255       sed_qc_const = k_st * ( 3.0_wp / ( 4.0_wp * pi * rho_l )                &
256                          )**( 2.0_wp / 3.0_wp ) *                             &
257                      EXP( 5.0_wp * LOG( sigma_gc )**2 )
258
259!
260!--    Calculate timestep according to precipitation
261       IF ( microphysics_seifert )  THEN
262          dt_precipitation = c_sedimentation * MINVAL( dzu(nzb+2:nzt) ) /      &
263                             w_precipitation
264       ENDIF
265
266!
267!--    Pre-calculate frequently calculated fractions of pi and rho_l
268       pirho_l  = pi * rho_l / 6.0_wp
269       dpirho_l = 1.0_wp / pirho_l
270
271!
272!--    Allocate 1D microphysics arrays
273       ALLOCATE ( nc_1d(nzb:nzt+1), pt_1d(nzb:nzt+1), q_1d(nzb:nzt+1),         &
274                  qc_1d(nzb:nzt+1) )
275
276       IF ( microphysics_seifert )  THEN
277          ALLOCATE ( nr_1d(nzb:nzt+1), qr_1d(nzb:nzt+1) )
278       ENDIF
279
280!
281!--    Initialze nc_1d with nc_const
282       nc_1d = nc_const
283
284    END SUBROUTINE microphysics_init
285
286
287!------------------------------------------------------------------------------!
288! Description:
289! ------------
290!> Control of microphysics for all grid points
291!------------------------------------------------------------------------------!
292    SUBROUTINE microphysics_control
293
294       USE arrays_3d,                                                          &
295           ONLY:  hyp, pt_init, prr, zu
296
297       USE cloud_parameters,                                                   &
298           ONLY:  cp, hyrho, pt_d_t, r_d, t_d_pt
299
300       USE control_parameters,                                                 &
301           ONLY:  call_microphysics_at_all_substeps, dt_3d, g,                 &
302                  intermediate_timestep_count, large_scale_forcing,            &
303                  lsf_surf, microphysics_kessler, microphysics_seifert,        &
304                  pt_surface, rho_surface,surface_pressure
305
306       USE indices,                                                            &
307           ONLY:  nzb, nzt
308
309       USE kinds
310
311       USE statistics,                                                         &
312           ONLY:  weight_pres
313
314       IMPLICIT NONE
315
316       INTEGER(iwp) ::  k                 !<
317
318       REAL(wp)     ::  t_surface         !<
319
320       IF ( large_scale_forcing  .AND.  lsf_surf ) THEN
321!
322!--       Calculate:
323!--       pt / t : ratio of potential and actual temperature (pt_d_t)
324!--       t / pt : ratio of actual and potential temperature (t_d_pt)
325!--       p_0(z) : vertical profile of the hydrostatic pressure (hyp)
326          t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**0.286_wp
327          DO  k = nzb, nzt+1
328             hyp(k)    = surface_pressure * 100.0_wp * &
329                         ( ( t_surface - g / cp * zu(k) ) /                    &
330                         t_surface )**(1.0_wp / 0.286_wp)
331             pt_d_t(k) = ( 100000.0_wp / hyp(k) )**0.286_wp
332             t_d_pt(k) = 1.0_wp / pt_d_t(k)
333             hyrho(k)  = hyp(k) / ( r_d * t_d_pt(k) * pt_init(k) )       
334          ENDDO
335
336!
337!--       Compute reference density
338          rho_surface = surface_pressure * 100.0_wp / ( r_d * t_surface )
339       ENDIF
340
341!
342!--    Compute length of time step
343       IF ( call_microphysics_at_all_substeps )  THEN
344          dt_micro = dt_3d * weight_pres(intermediate_timestep_count)
345       ELSE
346          dt_micro = dt_3d
347       ENDIF
348
349!
350!--    Reset precipitation rate
351       IF ( intermediate_timestep_count == 1 )  prr = 0.0_wp
352
353!
354!--    Compute cloud physics
355       IF ( microphysics_kessler )  THEN
356
357          CALL autoconversion_kessler
358          IF ( cloud_water_sedimentation )  CALL sedimentation_cloud
359
360       ELSEIF ( microphysics_seifert )  THEN
361
362          CALL adjust_cloud
363          CALL autoconversion
364          CALL accretion
365          CALL selfcollection_breakup
366          CALL evaporation_rain
367          CALL sedimentation_rain
368          IF ( cloud_water_sedimentation )  CALL sedimentation_cloud
369
370       ENDIF
371
372       CALL calc_precipitation_amount
373
374    END SUBROUTINE microphysics_control
375
376!------------------------------------------------------------------------------!
377! Description:
378! ------------
379!> Adjust number of raindrops to avoid nonlinear effects in sedimentation and
380!> evaporation of rain drops due to too small or too big weights
381!> of rain drops (Stevens and Seifert, 2008).
382!------------------------------------------------------------------------------!
383    SUBROUTINE adjust_cloud
384
385       USE arrays_3d,                                                          &
386           ONLY:  qr, nr
387
388       USE cloud_parameters,                                                   &
389           ONLY:  hyrho
390
391       USE cpulog,                                                             &
392           ONLY:  cpu_log, log_point_s
393
394       USE indices,                                                            &
395           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner, nzt
396
397       USE kinds
398
399       IMPLICIT NONE
400
401       INTEGER(iwp) ::  i                 !<
402       INTEGER(iwp) ::  j                 !<
403       INTEGER(iwp) ::  k                 !<
404
405       CALL cpu_log( log_point_s(54), 'adjust_cloud', 'start' )
406
407       DO  i = nxl, nxr
408          DO  j = nys, nyn
409             DO  k = nzb_s_inner(j,i)+1, nzt
410                IF ( qr(k,j,i) <= eps_sb )  THEN
411                   qr(k,j,i) = 0.0_wp
412                   nr(k,j,i) = 0.0_wp
413                ELSE
414                   IF ( nr(k,j,i) * xrmin > qr(k,j,i) * hyrho(k) )  THEN
415                      nr(k,j,i) = qr(k,j,i) * hyrho(k) / xrmin
416                   ELSEIF ( nr(k,j,i) * xrmax < qr(k,j,i) * hyrho(k) )  THEN
417                      nr(k,j,i) = qr(k,j,i) * hyrho(k) / xrmax
418                   ENDIF
419                ENDIF
420             ENDDO
421          ENDDO
422       ENDDO
423
424       CALL cpu_log( log_point_s(54), 'adjust_cloud', 'stop' )
425
426    END SUBROUTINE adjust_cloud
427
428
429!------------------------------------------------------------------------------!
430! Description:
431! ------------
432!> Autoconversion rate (Seifert and Beheng, 2006).
433!------------------------------------------------------------------------------!
434    SUBROUTINE autoconversion
435
436       USE arrays_3d,                                                          &
437           ONLY:  diss, dzu, nr, qc, qr
438
439       USE cloud_parameters,                                                   &
440           ONLY:  hyrho
441
442       USE control_parameters,                                                 &
443           ONLY:  rho_surface
444
445       USE cpulog,                                                             &
446           ONLY:  cpu_log, log_point_s
447
448       USE grid_variables,                                                     &
449           ONLY:  dx, dy
450
451       USE indices,                                                            &
452           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner, nzt
453
454       USE kinds
455
456       IMPLICIT NONE
457
458       INTEGER(iwp) ::  i                 !<
459       INTEGER(iwp) ::  j                 !<
460       INTEGER(iwp) ::  k                 !<
461
462       REAL(wp)     ::  alpha_cc          !<                   
463       REAL(wp)     ::  autocon           !<
464       REAL(wp)     ::  dissipation       !<
465       REAL(wp)     ::  k_au              !<
466       REAL(wp)     ::  l_mix             !<
467       REAL(wp)     ::  nu_c              !<
468       REAL(wp)     ::  phi_au            !<
469       REAL(wp)     ::  r_cc              !<
470       REAL(wp)     ::  rc                !<
471       REAL(wp)     ::  re_lambda         !<
472       REAL(wp)     ::  sigma_cc          !<
473       REAL(wp)     ::  tau_cloud         !<
474       REAL(wp)     ::  xc                !<
475
476       CALL cpu_log( log_point_s(55), 'autoconversion', 'start' )
477
478       DO  i = nxl, nxr
479          DO  j = nys, nyn
480             DO  k = nzb_s_inner(j,i)+1, nzt
481
482                IF ( qc(k,j,i) > eps_sb )  THEN
483
484                   k_au = k_cc / ( 20.0_wp * x0 )
485!
486!--                Intern time scale of coagulation (Seifert and Beheng, 2006):
487!--                (1.0_wp - qc(k,j,i) / ( qc(k,j,i) + qr(k,j,i) ))
488                   tau_cloud = 1.0_wp - qc(k,j,i) / ( qr(k,j,i) + qc(k,j,i) )
489!
490!--                Universal function for autoconversion process
491!--                (Seifert and Beheng, 2006):
492                   phi_au = 600.0_wp * tau_cloud**0.68_wp *                    &
493                            ( 1.0_wp - tau_cloud**0.68_wp )**3
494!
495!--                Shape parameter of gamma distribution (Geoffroy et al., 2010):
496!--                (Use constant nu_c = 1.0_wp instead?)
497                   nu_c = 1.0_wp !MAX( 0.0_wp, 1580.0_wp * hyrho(k) * qc(k,j,i) - 0.28_wp )
498!
499!--                Mean weight of cloud droplets:
500                   xc = hyrho(k) * qc(k,j,i) / nc_const
501!
502!--                Parameterized turbulence effects on autoconversion (Seifert,
503!--                Nuijens and Stevens, 2010)
504                   IF ( collision_turbulence )  THEN
505!
506!--                   Weight averaged radius of cloud droplets:
507                      rc = 0.5_wp * ( xc * dpirho_l )**( 1.0_wp / 3.0_wp )
508
509                      alpha_cc = ( a_1 + a_2 * nu_c ) / ( 1.0_wp + a_3 * nu_c )
510                      r_cc     = ( b_1 + b_2 * nu_c ) / ( 1.0_wp + b_3 * nu_c )
511                      sigma_cc = ( c_1 + c_2 * nu_c ) / ( 1.0_wp + c_3 * nu_c )
512!
513!--                   Mixing length (neglecting distance to ground and
514!--                   stratification)
515                      l_mix = ( dx * dy * dzu(k) )**( 1.0_wp / 3.0_wp )
516!
517!--                   Limit dissipation rate according to Seifert, Nuijens and
518!--                   Stevens (2010)
519                      dissipation = MIN( 0.06_wp, diss(k,j,i) )
520!
521!--                   Compute Taylor-microscale Reynolds number:
522                      re_lambda = 6.0_wp / 11.0_wp *                           &
523                                  ( l_mix / c_const )**( 2.0_wp / 3.0_wp ) *   &
524                                  SQRT( 15.0_wp / kin_vis_air ) *              &
525                                  dissipation**( 1.0_wp / 6.0_wp )
526!
527!--                   The factor of 1.0E4 is needed to convert the dissipation
528!--                   rate from m2 s-3 to cm2 s-3.
529                      k_au = k_au * ( 1.0_wp +                                 &
530                                      dissipation * 1.0E4_wp *                 &
531                                      ( re_lambda * 1.0E-3_wp )**0.25_wp *     &
532                                      ( alpha_cc * EXP( -1.0_wp * ( ( rc -     &     
533                                                                      r_cc ) / &
534                                                        sigma_cc )**2          &
535                                                      ) + beta_cc              &
536                                      )                                        &
537                                    )
538                   ENDIF
539!
540!--                Autoconversion rate (Seifert and Beheng, 2006):
541                   autocon = k_au * ( nu_c + 2.0_wp ) * ( nu_c + 4.0_wp ) /    &
542                             ( nu_c + 1.0_wp )**2 * qc(k,j,i)**2 * xc**2 *     &
543                             ( 1.0_wp + phi_au / ( 1.0_wp - tau_cloud )**2 ) * &
544                             rho_surface
545                   autocon = MIN( autocon, qc(k,j,i) / dt_micro )
546
547                   qr(k,j,i) = qr(k,j,i) + autocon * dt_micro
548                   qc(k,j,i) = qc(k,j,i) - autocon * dt_micro 
549                   nr(k,j,i) = nr(k,j,i) + autocon / x0 * hyrho(k) * dt_micro
550
551                ENDIF
552
553             ENDDO
554          ENDDO
555       ENDDO
556
557       CALL cpu_log( log_point_s(55), 'autoconversion', 'stop' )
558
559    END SUBROUTINE autoconversion
560
561
562!------------------------------------------------------------------------------!
563! Description:
564! ------------
565!> Autoconversion process (Kessler, 1969).
566!------------------------------------------------------------------------------!
567    SUBROUTINE autoconversion_kessler
568
569       USE arrays_3d,                                                          &
570           ONLY:  dzw, pt, prr, q, qc
571
572       USE cloud_parameters,                                                   &
573           ONLY:  l_d_cp, pt_d_t
574
575       USE indices,                                                            &
576           ONLY:  nxl, nxr, nyn, nys, nzb_s_inner, nzt
577
578       USE kinds
579
580
581       IMPLICIT NONE
582
583       INTEGER(iwp) ::  i !<
584       INTEGER(iwp) ::  j !<
585       INTEGER(iwp) ::  k !<
586
587       REAL(wp)    ::  dqdt_precip !<
588
589       DO  i = nxl, nxr
590          DO  j = nys, nyn
591             DO  k = nzb_s_inner(j,i)+1, nzt
592
593                IF ( qc(k,j,i) > ql_crit )  THEN
594                   dqdt_precip = prec_time_const * ( qc(k,j,i) - ql_crit )
595                ELSE
596                   dqdt_precip = 0.0_wp
597                ENDIF
598
599                qc(k,j,i) = qc(k,j,i) - dqdt_precip * dt_micro
600                q(k,j,i)  = q(k,j,i)  - dqdt_precip * dt_micro
601                pt(k,j,i) = pt(k,j,i) + dqdt_precip * dt_micro * l_d_cp *      &
602                                        pt_d_t(k)
603
604!
605!--             Compute the rain rate (stored on surface grid point)
606                prr(nzb_s_inner(j,i),j,i) = prr(nzb_s_inner(j,i),j,i) +        &
607                                            dqdt_precip * dzw(k)
608
609             ENDDO
610          ENDDO
611       ENDDO
612
613   END SUBROUTINE autoconversion_kessler
614
615
616!------------------------------------------------------------------------------!
617! Description:
618! ------------
619!> Accretion rate (Seifert and Beheng, 2006).
620!------------------------------------------------------------------------------!
621    SUBROUTINE accretion
622
623       USE arrays_3d,                                                          &
624           ONLY:  diss, qc, qr
625
626       USE cloud_parameters,                                                   &
627           ONLY:  hyrho
628
629       USE control_parameters,                                                 &
630           ONLY:  rho_surface
631
632       USE cpulog,                                                             &
633           ONLY:  cpu_log, log_point_s
634
635       USE indices,                                                            &
636           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner, nzt
637
638       USE kinds
639
640       IMPLICIT NONE
641
642       INTEGER(iwp) ::  i                 !<
643       INTEGER(iwp) ::  j                 !<
644       INTEGER(iwp) ::  k                 !<
645
646       REAL(wp)     ::  accr              !<
647       REAL(wp)     ::  k_cr              !<
648       REAL(wp)     ::  phi_ac            !<
649       REAL(wp)     ::  tau_cloud         !<
650
651       CALL cpu_log( log_point_s(56), 'accretion', 'start' )
652
653       DO  i = nxl, nxr
654          DO  j = nys, nyn
655             DO  k = nzb_s_inner(j,i)+1, nzt
656
657                IF ( ( qc(k,j,i) > eps_sb )  .AND.  ( qr(k,j,i) > eps_sb ) )  THEN
658!
659!--                Intern time scale of coagulation (Seifert and Beheng, 2006):
660                   tau_cloud = 1.0_wp - qc(k,j,i) / ( qc(k,j,i) + qr(k,j,i) ) 
661!
662!--                Universal function for accretion process (Seifert and
663!--                Beheng, 2001):
664                   phi_ac = ( tau_cloud / ( tau_cloud + 5.0E-5_wp ) )**4
665!
666!--                Parameterized turbulence effects on autoconversion (Seifert,
667!--                Nuijens and Stevens, 2010). The factor of 1.0E4 is needed to
668!--                convert the dissipation rate (diss) from m2 s-3 to cm2 s-3.
669                   IF ( collision_turbulence )  THEN
670                      k_cr = k_cr0 * ( 1.0_wp + 0.05_wp *                      &
671                                       MIN( 600.0_wp,                          &
672                                            diss(k,j,i) * 1.0E4_wp )**0.25_wp  &
673                                     )
674                   ELSE
675                      k_cr = k_cr0                       
676                   ENDIF
677!
678!--                Accretion rate (Seifert and Beheng, 2006):
679                   accr = k_cr * qc(k,j,i) * qr(k,j,i) * phi_ac *              &
680                          SQRT( rho_surface * hyrho(k) )
681                   accr = MIN( accr, qc(k,j,i) / dt_micro )
682
683                   qr(k,j,i) = qr(k,j,i) + accr * dt_micro 
684                   qc(k,j,i) = qc(k,j,i) - accr * dt_micro
685
686                ENDIF
687
688             ENDDO
689          ENDDO
690       ENDDO
691
692       CALL cpu_log( log_point_s(56), 'accretion', 'stop' )
693
694    END SUBROUTINE accretion
695
696
697!------------------------------------------------------------------------------!
698! Description:
699! ------------
700!> Collisional breakup rate (Seifert, 2008).
701!------------------------------------------------------------------------------!
702    SUBROUTINE selfcollection_breakup
703
704       USE arrays_3d,                                                          &
705           ONLY:  nr, qr
706
707       USE cloud_parameters,                                                   &
708           ONLY:  hyrho
709
710       USE control_parameters,                                                 &
711           ONLY:  rho_surface
712
713       USE cpulog,                                                             &
714           ONLY:  cpu_log, log_point_s
715
716       USE indices,                                                            &
717           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner, nzt
718
719       USE kinds
720   
721       IMPLICIT NONE
722
723       INTEGER(iwp) ::  i                 !<
724       INTEGER(iwp) ::  j                 !<
725       INTEGER(iwp) ::  k                 !<
726
727       REAL(wp)     ::  breakup           !<
728       REAL(wp)     ::  dr                !<
729       REAL(wp)     ::  phi_br            !<
730       REAL(wp)     ::  selfcoll          !<
731
732       CALL cpu_log( log_point_s(57), 'selfcollection', 'start' )
733
734       DO  i = nxl, nxr
735          DO  j = nys, nyn
736             DO  k = nzb_s_inner(j,i)+1, nzt
737                IF ( qr(k,j,i) > eps_sb )  THEN
738!
739!--                Selfcollection rate (Seifert and Beheng, 2001):
740                   selfcoll = k_rr * nr(k,j,i) * qr(k,j,i) *                   &
741                              SQRT( hyrho(k) * rho_surface )
742!
743!--                Weight averaged diameter of rain drops:
744                   dr = ( hyrho(k) * qr(k,j,i) /                               &
745                          nr(k,j,i) * dpirho_l )**( 1.0_wp / 3.0_wp )
746!
747!--                Collisional breakup rate (Seifert, 2008):
748                   IF ( dr >= 0.3E-3_wp )  THEN
749                      phi_br  = k_br * ( dr - 1.1E-3_wp )
750                      breakup = selfcoll * ( phi_br + 1.0_wp )
751                   ELSE
752                      breakup = 0.0_wp
753                   ENDIF
754
755                   selfcoll = MAX( breakup - selfcoll, -nr(k,j,i) / dt_micro )
756                   nr(k,j,i) = nr(k,j,i) + selfcoll * dt_micro
757
758                ENDIF         
759             ENDDO
760          ENDDO
761       ENDDO
762
763       CALL cpu_log( log_point_s(57), 'selfcollection', 'stop' )
764
765    END SUBROUTINE selfcollection_breakup
766
767
768!------------------------------------------------------------------------------!
769! Description:
770! ------------
771!> Evaporation of precipitable water. Condensation is neglected for
772!> precipitable water.
773!------------------------------------------------------------------------------!
774    SUBROUTINE evaporation_rain
775
776       USE arrays_3d,                                                          &
777           ONLY:  hyp, nr, pt, q,  qc, qr
778
779       USE cloud_parameters,                                                   &
780           ONLY:  hyrho, l_d_cp, l_d_r, l_v, r_v, t_d_pt
781
782       USE constants,                                                          &
783           ONLY:  pi
784
785       USE cpulog,                                                             &
786           ONLY:  cpu_log, log_point_s
787
788       USE indices,                                                            &
789           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner, nzt
790
791       USE kinds
792
793       IMPLICIT NONE
794
795       INTEGER(iwp) ::  i                 !<
796       INTEGER(iwp) ::  j                 !<
797       INTEGER(iwp) ::  k                 !<
798
799       REAL(wp)     ::  alpha             !<
800       REAL(wp)     ::  dr                !<
801       REAL(wp)     ::  e_s               !<
802       REAL(wp)     ::  evap              !<
803       REAL(wp)     ::  evap_nr           !<
804       REAL(wp)     ::  f_vent            !<
805       REAL(wp)     ::  g_evap            !<
806       REAL(wp)     ::  lambda_r          !<
807       REAL(wp)     ::  mu_r              !<
808       REAL(wp)     ::  mu_r_2            !<
809       REAL(wp)     ::  mu_r_5d2          !<
810       REAL(wp)     ::  nr_0              !<
811       REAL(wp)     ::  q_s               !<
812       REAL(wp)     ::  sat               !<
813       REAL(wp)     ::  t_l               !<
814       REAL(wp)     ::  temp              !<
815       REAL(wp)     ::  xr                !<
816
817       CALL cpu_log( log_point_s(58), 'evaporation', 'start' )
818
819       DO  i = nxl, nxr
820          DO  j = nys, nyn
821             DO  k = nzb_s_inner(j,i)+1, nzt
822                IF ( qr(k,j,i) > eps_sb )  THEN
823!
824!--                Actual liquid water temperature:
825                   t_l = t_d_pt(k) * pt(k,j,i)
826!
827!--                Saturation vapor pressure at t_l:
828                   e_s = 610.78_wp * EXP( 17.269_wp * ( t_l - 273.16_wp ) /    &
829                                          ( t_l - 35.86_wp )                   &
830                                        )
831!
832!--                Computation of saturation humidity:
833                   q_s   = 0.622_wp * e_s / ( hyp(k) - 0.378_wp * e_s )
834                   alpha = 0.622_wp * l_d_r * l_d_cp / ( t_l * t_l )
835                   q_s   = q_s * ( 1.0_wp + alpha * q(k,j,i) ) /               &
836                           ( 1.0_wp + alpha * q_s )
837!
838!--                Supersaturation:
839                   sat   = ( q(k,j,i) - qr(k,j,i) - qc(k,j,i) ) / q_s - 1.0_wp
840!
841!--                Evaporation needs only to be calculated in subsaturated regions
842                   IF ( sat < 0.0_wp )  THEN
843!
844!--                   Actual temperature:
845                      temp = t_l + l_d_cp * ( qc(k,j,i) + qr(k,j,i) )
846             
847                      g_evap = 1.0_wp / ( ( l_v / ( r_v * temp ) - 1.0_wp ) *     &
848                                          l_v / ( thermal_conductivity_l * temp ) &
849                                          + r_v * temp / ( diff_coeff_l * e_s )   &
850                                        )
851!
852!--                   Mean weight of rain drops
853                      xr = hyrho(k) * qr(k,j,i) / nr(k,j,i)
854!
855!--                   Weight averaged diameter of rain drops:
856                      dr = ( xr * dpirho_l )**( 1.0_wp / 3.0_wp )
857!
858!--                   Compute ventilation factor and intercept parameter
859!--                   (Seifert and Beheng, 2006; Seifert, 2008):
860                      IF ( ventilation_effect )  THEN
861!
862!--                      Shape parameter of gamma distribution (Milbrandt and Yau,
863!--                      2005; Stevens and Seifert, 2008):
864                         mu_r = 10.0_wp * ( 1.0_wp + TANH( 1.2E3_wp *          &
865                                                          ( dr - 1.4E-3_wp ) ) )
866!
867!--                      Slope parameter of gamma distribution (Seifert, 2008):
868                         lambda_r = ( ( mu_r + 3.0_wp ) * ( mu_r + 2.0_wp ) *  &
869                                      ( mu_r + 1.0_wp )                        &
870                                    )**( 1.0_wp / 3.0_wp ) / dr
871
872                         mu_r_2   = mu_r + 2.0_wp
873                         mu_r_5d2 = mu_r + 2.5_wp 
874
875                         f_vent = a_vent * gamm( mu_r_2 ) *                     &
876                                  lambda_r**( -mu_r_2 ) + b_vent *              &
877                                  schmidt_p_1d3 * SQRT( a_term / kin_vis_air ) *&
878                                  gamm( mu_r_5d2 ) * lambda_r**( -mu_r_5d2 ) *  &
879                                  ( 1.0_wp -                                    &
880                                    0.5_wp * ( b_term / a_term ) *              &
881                                    ( lambda_r / ( c_term + lambda_r )          &
882                                    )**mu_r_5d2 -                               &
883                                    0.125_wp * ( b_term / a_term )**2 *         &
884                                    ( lambda_r / ( 2.0_wp * c_term + lambda_r ) &
885                                    )**mu_r_5d2 -                               &
886                                    0.0625_wp * ( b_term / a_term )**3 *        &
887                                    ( lambda_r / ( 3.0_wp * c_term + lambda_r ) &
888                                    )**mu_r_5d2 -                               &
889                                    0.0390625_wp * ( b_term / a_term )**4 *     & 
890                                    ( lambda_r / ( 4.0_wp * c_term + lambda_r ) &
891                                    )**mu_r_5d2                                 &
892                                  )
893
894                         nr_0   = nr(k,j,i) * lambda_r**( mu_r + 1.0_wp ) /    &
895                                  gamm( mu_r + 1.0_wp ) 
896                      ELSE
897                         f_vent = 1.0_wp
898                         nr_0   = nr(k,j,i) * dr
899                      ENDIF
900   !
901   !--                Evaporation rate of rain water content (Seifert and
902   !--                Beheng, 2006):
903                      evap    = 2.0_wp * pi * nr_0 * g_evap * f_vent * sat /   &
904                                hyrho(k)
905                      evap    = MAX( evap, -qr(k,j,i) / dt_micro )
906                      evap_nr = MAX( c_evap * evap / xr * hyrho(k),            &
907                                     -nr(k,j,i) / dt_micro )
908
909                      qr(k,j,i) = qr(k,j,i) + evap    * dt_micro
910                      nr(k,j,i) = nr(k,j,i) + evap_nr * dt_micro
911
912                   ENDIF
913                ENDIF         
914
915             ENDDO
916          ENDDO
917       ENDDO
918
919       CALL cpu_log( log_point_s(58), 'evaporation', 'stop' )
920
921    END SUBROUTINE evaporation_rain
922
923
924!------------------------------------------------------------------------------!
925! Description:
926! ------------
927!> Sedimentation of cloud droplets (Ackermann et al., 2009, MWR).
928!------------------------------------------------------------------------------!
929    SUBROUTINE sedimentation_cloud
930
931       USE arrays_3d,                                                          &
932           ONLY:  ddzu, dzu, pt, prr, q, qc
933
934       USE cloud_parameters,                                                   &
935           ONLY:  hyrho, l_d_cp, pt_d_t
936
937       USE control_parameters,                                                 &
938           ONLY:  call_microphysics_at_all_substeps, intermediate_timestep_count
939
940       USE cpulog,                                                             &
941           ONLY:  cpu_log, log_point_s
942
943       USE indices,                                                            &
944           ONLY:  nxl, nxr, nys, nyn, nzb, nzb_s_inner, nzt
945
946       USE kinds
947
948       USE statistics,                                                         &
949           ONLY:  weight_substep
950   
951
952       IMPLICIT NONE
953
954       INTEGER(iwp) ::  i             !<
955       INTEGER(iwp) ::  j             !<
956       INTEGER(iwp) ::  k             !<
957
958       REAL(wp), DIMENSION(nzb:nzt+1) :: sed_qc !<
959
960       CALL cpu_log( log_point_s(59), 'sed_cloud', 'start' )
961
962       sed_qc(nzt+1) = 0.0_wp
963
964       DO  i = nxl, nxr
965          DO  j = nys, nyn
966             DO  k = nzt, nzb_s_inner(j,i)+1, -1
967
968                IF ( qc(k,j,i) > eps_sb )  THEN
969                   sed_qc(k) = sed_qc_const * nc_const**( -2.0_wp / 3.0_wp ) * &
970                               ( qc(k,j,i) * hyrho(k) )**( 5.0_wp / 3.0_wp )
971                ELSE
972                   sed_qc(k) = 0.0_wp
973                ENDIF
974
975                sed_qc(k) = MIN( sed_qc(k), hyrho(k) * dzu(k+1) * q(k,j,i) /   &
976                                            dt_micro + sed_qc(k+1)             &
977                               )
978
979                q(k,j,i)  = q(k,j,i)  + ( sed_qc(k+1) - sed_qc(k) ) *          &
980                                        ddzu(k+1) / hyrho(k) * dt_micro
981                qc(k,j,i) = qc(k,j,i) + ( sed_qc(k+1) - sed_qc(k) ) *          &
982                                        ddzu(k+1) / hyrho(k) * dt_micro
983                pt(k,j,i) = pt(k,j,i) - ( sed_qc(k+1) - sed_qc(k) ) *          &
984                                        ddzu(k+1) / hyrho(k) * l_d_cp *        &
985                                        pt_d_t(k) * dt_micro
986
987!
988!--             Compute the precipitation rate due to cloud (fog) droplets
989                IF ( call_microphysics_at_all_substeps )  THEN
990                   prr(k,j,i) = prr(k,j,i) +  sed_qc(k) / hyrho(k)             &
991                                * weight_substep(intermediate_timestep_count)
992                ELSE
993                   prr(k,j,i) = prr(k,j,i) + sed_qc(k) / hyrho(k)
994                ENDIF
995
996             ENDDO
997          ENDDO
998       ENDDO
999
1000       CALL cpu_log( log_point_s(59), 'sed_cloud', 'stop' )
1001
1002    END SUBROUTINE sedimentation_cloud
1003
1004
1005!------------------------------------------------------------------------------!
1006! Description:
1007! ------------
1008!> Computation of sedimentation flux. Implementation according to Stevens
1009!> and Seifert (2008). Code is based on UCLA-LES.
1010!------------------------------------------------------------------------------!
1011    SUBROUTINE sedimentation_rain
1012
1013       USE arrays_3d,                                                          &
1014           ONLY:  ddzu, dzu, nr, pt, prr, q, qr
1015
1016       USE cloud_parameters,                                                   &
1017           ONLY:  hyrho, l_d_cp, pt_d_t
1018
1019       USE control_parameters,                                                 &
1020           ONLY:  call_microphysics_at_all_substeps, intermediate_timestep_count
1021       USE cpulog,                                                             &
1022           ONLY:  cpu_log, log_point_s
1023
1024       USE indices,                                                            &
1025           ONLY:  nxl, nxr, nys, nyn, nzb, nzb_s_inner, nzt
1026
1027       USE kinds
1028
1029       USE statistics,                                                         &
1030           ONLY:  weight_substep
1031       
1032       IMPLICIT NONE
1033
1034       INTEGER(iwp) ::  i                          !<
1035       INTEGER(iwp) ::  j                          !<
1036       INTEGER(iwp) ::  k                          !<
1037       INTEGER(iwp) ::  k_run                      !<
1038
1039       REAL(wp)     ::  c_run                      !<
1040       REAL(wp)     ::  d_max                      !<
1041       REAL(wp)     ::  d_mean                     !<
1042       REAL(wp)     ::  d_min                      !<
1043       REAL(wp)     ::  dr                         !<
1044       REAL(wp)     ::  flux                       !<
1045       REAL(wp)     ::  lambda_r                   !<
1046       REAL(wp)     ::  mu_r                       !<
1047       REAL(wp)     ::  z_run                      !<
1048
1049       REAL(wp), DIMENSION(nzb:nzt+1) ::  c_nr     !<
1050       REAL(wp), DIMENSION(nzb:nzt+1) ::  c_qr     !<
1051       REAL(wp), DIMENSION(nzb:nzt+1) ::  nr_slope !<
1052       REAL(wp), DIMENSION(nzb:nzt+1) ::  qr_slope !<
1053       REAL(wp), DIMENSION(nzb:nzt+1) ::  sed_nr   !<
1054       REAL(wp), DIMENSION(nzb:nzt+1) ::  sed_qr   !<
1055       REAL(wp), DIMENSION(nzb:nzt+1) ::  w_nr     !<
1056       REAL(wp), DIMENSION(nzb:nzt+1) ::  w_qr     !<
1057
1058       CALL cpu_log( log_point_s(60), 'sed_rain', 'start' )
1059
1060!
1061!--    Compute velocities
1062       DO  i = nxl, nxr
1063          DO  j = nys, nyn
1064             DO  k = nzb_s_inner(j,i)+1, nzt
1065                IF ( qr(k,j,i) > eps_sb )  THEN
1066!
1067!--                Weight averaged diameter of rain drops:
1068                   dr = ( hyrho(k) * qr(k,j,i) /                               &
1069                          nr(k,j,i) * dpirho_l )**( 1.0_wp / 3.0_wp )
1070!
1071!--                Shape parameter of gamma distribution (Milbrandt and Yau, 2005;
1072!--                Stevens and Seifert, 2008):
1073                   mu_r = 10.0_wp * ( 1.0_wp + TANH( 1.2E3_wp *                &
1074                                                     ( dr - 1.4E-3_wp ) ) )
1075!
1076!--                Slope parameter of gamma distribution (Seifert, 2008):
1077                   lambda_r = ( ( mu_r + 3.0_wp ) * ( mu_r + 2.0_wp ) *        &
1078                                ( mu_r + 1.0_wp ) )**( 1.0_wp / 3.0_wp ) / dr
1079
1080                   w_nr(k) = MAX( 0.1_wp, MIN( 20.0_wp,                        &
1081                                               a_term - b_term * ( 1.0_wp +    &
1082                                                  c_term /                     &
1083                                                  lambda_r )**( -1.0_wp *      &
1084                                                  ( mu_r + 1.0_wp ) )          &
1085                                              )                                &
1086                                )
1087
1088                   w_qr(k) = MAX( 0.1_wp, MIN( 20.0_wp,                        &
1089                                               a_term - b_term * ( 1.0_wp +    &
1090                                                  c_term /                     &
1091                                                  lambda_r )**( -1.0_wp *      &
1092                                                  ( mu_r + 4.0_wp ) )          &
1093                                             )                                 &
1094                                )
1095                ELSE
1096                   w_nr(k) = 0.0_wp
1097                   w_qr(k) = 0.0_wp
1098                ENDIF
1099             ENDDO
1100!
1101!--          Adjust boundary values
1102             w_nr(nzb_s_inner(j,i)) = w_nr(nzb_s_inner(j,i)+1)
1103             w_qr(nzb_s_inner(j,i)) = w_qr(nzb_s_inner(j,i)+1)
1104             w_nr(nzt+1) = 0.0_wp
1105             w_qr(nzt+1) = 0.0_wp
1106!
1107!--          Compute Courant number
1108             DO  k = nzb_s_inner(j,i)+1, nzt
1109                c_nr(k) = 0.25_wp * ( w_nr(k-1) +                              &
1110                                      2.0_wp * w_nr(k) + w_nr(k+1) ) *         &
1111                          dt_micro * ddzu(k)
1112                c_qr(k) = 0.25_wp * ( w_qr(k-1) +                              &
1113                                      2.0_wp * w_qr(k) + w_qr(k+1) ) *         &
1114                          dt_micro * ddzu(k)
1115             ENDDO     
1116!
1117!--          Limit slopes with monotonized centered (MC) limiter (van Leer, 1977):
1118             IF ( limiter_sedimentation )  THEN
1119
1120                DO k = nzb_s_inner(j,i)+1, nzt
1121                   d_mean = 0.5_wp * ( qr(k+1,j,i) - qr(k-1,j,i) )
1122                   d_min  = qr(k,j,i) - MIN( qr(k+1,j,i), qr(k,j,i), qr(k-1,j,i) )
1123                   d_max  = MAX( qr(k+1,j,i), qr(k,j,i), qr(k-1,j,i) ) - qr(k,j,i)
1124
1125                   qr_slope(k) = SIGN(1.0_wp, d_mean) * MIN ( 2.0_wp * d_min,  &
1126                                                              2.0_wp * d_max,  &
1127                                                              ABS( d_mean ) )
1128
1129                   d_mean = 0.5_wp * ( nr(k+1,j,i) - nr(k-1,j,i) )
1130                   d_min  = nr(k,j,i) - MIN( nr(k+1,j,i), nr(k,j,i), nr(k-1,j,i) )
1131                   d_max  = MAX( nr(k+1,j,i), nr(k,j,i), nr(k-1,j,i) ) - nr(k,j,i)
1132
1133                   nr_slope(k) = SIGN(1.0_wp, d_mean) * MIN ( 2.0_wp * d_min,  &
1134                                                              2.0_wp * d_max,  &
1135                                                              ABS( d_mean ) )
1136                ENDDO
1137
1138             ELSE
1139
1140                nr_slope = 0.0_wp
1141                qr_slope = 0.0_wp
1142
1143             ENDIF
1144
1145             sed_nr(nzt+1) = 0.0_wp
1146             sed_qr(nzt+1) = 0.0_wp
1147!
1148!--          Compute sedimentation flux
1149             DO  k = nzt, nzb_s_inner(j,i)+1, -1
1150!
1151!--             Sum up all rain drop number densities which contribute to the flux
1152!--             through k-1/2
1153                flux  = 0.0_wp
1154                z_run = 0.0_wp ! height above z(k)
1155                k_run = k
1156                c_run = MIN( 1.0_wp, c_nr(k) )
1157                DO WHILE ( c_run > 0.0_wp  .AND.  k_run <= nzt )
1158                   flux  = flux + hyrho(k_run) *                               &
1159                           ( nr(k_run,j,i) + nr_slope(k_run) *                 &
1160                           ( 1.0_wp - c_run ) * 0.5_wp ) * c_run * dzu(k_run)
1161                   z_run = z_run + dzu(k_run)
1162                   k_run = k_run + 1
1163                   c_run = MIN( 1.0_wp, c_nr(k_run) - z_run * ddzu(k_run) )
1164                ENDDO
1165!
1166!--             It is not allowed to sediment more rain drop number density than
1167!--             available
1168                flux = MIN( flux,                                              &
1169                            hyrho(k) * dzu(k+1) * nr(k,j,i) + sed_nr(k+1) *    &
1170                            dt_micro                                           &
1171                          )
1172
1173                sed_nr(k) = flux / dt_micro
1174                nr(k,j,i) = nr(k,j,i) + ( sed_nr(k+1) - sed_nr(k) ) *          &
1175                                        ddzu(k+1) / hyrho(k) * dt_micro
1176!
1177!--             Sum up all rain water content which contributes to the flux
1178!--             through k-1/2
1179                flux  = 0.0_wp
1180                z_run = 0.0_wp ! height above z(k)
1181                k_run = k
1182                c_run = MIN( 1.0_wp, c_qr(k) )
1183
1184                DO WHILE ( c_run > 0.0_wp  .AND.  k_run <= nzt )
1185
1186                   flux  = flux + hyrho(k_run) * ( qr(k_run,j,i) +             &
1187                                  qr_slope(k_run) * ( 1.0_wp - c_run ) *       &
1188                                  0.5_wp ) * c_run * dzu(k_run)
1189                   z_run = z_run + dzu(k_run)
1190                   k_run = k_run + 1
1191                   c_run = MIN( 1.0_wp, c_qr(k_run) - z_run * ddzu(k_run) )
1192
1193                ENDDO
1194!
1195!--             It is not allowed to sediment more rain water content than
1196!--             available
1197                flux = MIN( flux,                                              &
1198                            hyrho(k) * dzu(k) * qr(k,j,i) + sed_qr(k+1) *      &
1199                            dt_micro                                           &
1200                          )
1201
1202                sed_qr(k) = flux / dt_micro
1203
1204                qr(k,j,i) = qr(k,j,i) + ( sed_qr(k+1) - sed_qr(k) ) *          &
1205                                        ddzu(k+1) / hyrho(k) * dt_micro
1206                q(k,j,i)  = q(k,j,i)  + ( sed_qr(k+1) - sed_qr(k) ) *          &
1207                                        ddzu(k+1) / hyrho(k) * dt_micro 
1208                pt(k,j,i) = pt(k,j,i) - ( sed_qr(k+1) - sed_qr(k) ) *          &
1209                                        ddzu(k+1) / hyrho(k) * l_d_cp *        &
1210                                        pt_d_t(k) * dt_micro
1211!
1212!--             Compute the rain rate
1213                IF ( call_microphysics_at_all_substeps )  THEN
1214                   prr(k,j,i) = prr(k,j,i) + sed_qr(k) / hyrho(k)             &
1215                                * weight_substep(intermediate_timestep_count)
1216                ELSE
1217                   prr(k,j,i) = prr(k,j,i) + sed_qr(k) / hyrho(k)
1218                ENDIF
1219
1220             ENDDO
1221          ENDDO
1222       ENDDO
1223
1224       CALL cpu_log( log_point_s(60), 'sed_rain', 'stop' )
1225
1226    END SUBROUTINE sedimentation_rain
1227
1228
1229!------------------------------------------------------------------------------!
1230! Description:
1231! ------------
1232!> Computation of the precipitation amount due to gravitational settling of
1233!> rain and cloud (fog) droplets
1234!------------------------------------------------------------------------------!
1235    SUBROUTINE calc_precipitation_amount
1236
1237       USE arrays_3d,                                                          &
1238           ONLY:  precipitation_amount, prr
1239
1240       USE cloud_parameters,                                                   &
1241           ONLY:  hyrho
1242
1243       USE control_parameters,                                                 &
1244           ONLY:  call_microphysics_at_all_substeps, dt_do2d_xy, dt_3d,        &
1245                  intermediate_timestep_count, intermediate_timestep_count_max,&
1246                  precipitation_amount_interval, time_do2d_xy
1247
1248       USE indices,                                                            &
1249           ONLY:  nxl, nxr, nys, nyn, nzb_s_inner
1250
1251       USE kinds
1252
1253       IMPLICIT NONE
1254
1255       INTEGER(iwp) ::  i                          !:
1256       INTEGER(iwp) ::  j                          !:
1257
1258
1259       IF ( ( dt_do2d_xy - time_do2d_xy ) < precipitation_amount_interval .AND.&
1260            ( .NOT. call_microphysics_at_all_substeps .OR.                     &
1261            intermediate_timestep_count == intermediate_timestep_count_max ) ) &
1262       THEN
1263
1264          DO  i = nxl, nxr
1265             DO  j = nys, nyn
1266
1267                precipitation_amount(j,i) = precipitation_amount(j,i) +        &
1268                                            prr(nzb_s_inner(j,i)+1,j,i) *      &
1269                                            hyrho(nzb_s_inner(j,i)+1) * dt_3d
1270
1271             ENDDO
1272          ENDDO
1273       ENDIF
1274
1275    END SUBROUTINE calc_precipitation_amount
1276
1277
1278!------------------------------------------------------------------------------!
1279! Description:
1280! ------------
1281!> Control of microphysics for grid points i,j
1282!------------------------------------------------------------------------------!
1283
1284    SUBROUTINE microphysics_control_ij( i, j )
1285
1286       USE arrays_3d,                                                          &
1287           ONLY:  hyp, nr, pt, pt_init, prr, q, qc, qr, zu
1288
1289       USE cloud_parameters,                                                   &
1290           ONLY:  cp, hyrho, pt_d_t, r_d, t_d_pt
1291
1292       USE control_parameters,                                                 &
1293           ONLY:  call_microphysics_at_all_substeps, dt_3d, g,                 &
1294                  intermediate_timestep_count, large_scale_forcing,            &
1295                  lsf_surf, microphysics_seifert, microphysics_kessler,        &
1296                  pt_surface, rho_surface, surface_pressure
1297
1298       USE indices,                                                            &
1299           ONLY:  nzb, nzt
1300
1301       USE kinds
1302
1303       USE statistics,                                                         &
1304           ONLY:  weight_pres
1305
1306       IMPLICIT NONE
1307
1308       INTEGER(iwp) ::  i                 !<
1309       INTEGER(iwp) ::  j                 !<
1310       INTEGER(iwp) ::  k                 !<
1311
1312       REAL(wp)     ::  t_surface         !<
1313
1314       IF ( large_scale_forcing  .AND.  lsf_surf ) THEN
1315!
1316!--       Calculate:
1317!--       pt / t : ratio of potential and actual temperature (pt_d_t)
1318!--       t / pt : ratio of actual and potential temperature (t_d_pt)
1319!--       p_0(z) : vertical profile of the hydrostatic pressure (hyp)
1320          t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**0.286_wp
1321          DO  k = nzb, nzt+1
1322             hyp(k)    = surface_pressure * 100.0_wp * &
1323                         ( ( t_surface - g / cp * zu(k) ) / t_surface )**(1.0_wp / 0.286_wp)
1324             pt_d_t(k) = ( 100000.0_wp / hyp(k) )**0.286_wp
1325             t_d_pt(k) = 1.0_wp / pt_d_t(k)
1326             hyrho(k)  = hyp(k) / ( r_d * t_d_pt(k) * pt_init(k) )       
1327          ENDDO
1328!
1329!--       Compute reference density
1330          rho_surface = surface_pressure * 100.0_wp / ( r_d * t_surface )
1331       ENDIF
1332
1333!
1334!--    Compute length of time step
1335       IF ( call_microphysics_at_all_substeps )  THEN
1336          dt_micro = dt_3d * weight_pres(intermediate_timestep_count)
1337       ELSE
1338          dt_micro = dt_3d
1339       ENDIF
1340
1341!
1342!--    Use 1d arrays
1343       q_1d(:)  = q(:,j,i)
1344       pt_1d(:) = pt(:,j,i)
1345       qc_1d(:) = qc(:,j,i)
1346       nc_1d(:) = nc_const
1347       IF ( microphysics_seifert )  THEN
1348          qr_1d(:) = qr(:,j,i)
1349          nr_1d(:) = nr(:,j,i)
1350       ENDIF
1351
1352!
1353!--    Reset precipitation rate
1354       IF ( intermediate_timestep_count == 1 )  prr(:,j,i) = 0.0_wp
1355
1356!
1357!--    Compute cloud physics
1358       IF( microphysics_kessler )  THEN
1359
1360          CALL autoconversion_kessler( i,j )
1361          IF ( cloud_water_sedimentation )  CALL sedimentation_cloud( i,j )
1362
1363       ELSEIF ( microphysics_seifert )  THEN
1364
1365          CALL adjust_cloud( i,j )
1366          CALL autoconversion( i,j )
1367          CALL accretion( i,j )
1368          CALL selfcollection_breakup( i,j )
1369          CALL evaporation_rain( i,j )
1370          CALL sedimentation_rain( i,j )
1371          IF ( cloud_water_sedimentation )  CALL sedimentation_cloud( i,j )
1372
1373       ENDIF
1374
1375       CALL calc_precipitation_amount( i,j )
1376
1377!
1378!--    Store results on the 3d arrays
1379       q(:,j,i)  = q_1d(:)
1380       pt(:,j,i) = pt_1d(:)
1381       IF ( microphysics_seifert )  THEN
1382          qr(:,j,i) = qr_1d(:)
1383          nr(:,j,i) = nr_1d(:)
1384       ENDIF
1385
1386    END SUBROUTINE microphysics_control_ij
1387
1388!------------------------------------------------------------------------------!
1389! Description:
1390! ------------
1391!> Adjust number of raindrops to avoid nonlinear effects in
1392!> sedimentation and evaporation of rain drops due to too small or
1393!> too big weights of rain drops (Stevens and Seifert, 2008).
1394!> The same procedure is applied to cloud droplets if they are determined
1395!> prognostically. Call for grid point i,j
1396!------------------------------------------------------------------------------!
1397    SUBROUTINE adjust_cloud_ij( i, j )
1398
1399       USE cloud_parameters,                                                   &
1400           ONLY:  hyrho
1401
1402       USE indices,                                                            &
1403           ONLY:  nzb_s_inner, nzt
1404
1405       USE kinds
1406
1407       IMPLICIT NONE
1408
1409       INTEGER(iwp) ::  i                 !<
1410       INTEGER(iwp) ::  j                 !<
1411       INTEGER(iwp) ::  k                 !<
1412
1413       DO  k = nzb_s_inner(j,i)+1, nzt
1414
1415          IF ( qr_1d(k) <= eps_sb )  THEN
1416             qr_1d(k) = 0.0_wp
1417             nr_1d(k) = 0.0_wp
1418          ELSE
1419!
1420!--          Adjust number of raindrops to avoid nonlinear effects in
1421!--          sedimentation and evaporation of rain drops due to too small or
1422!--          too big weights of rain drops (Stevens and Seifert, 2008).
1423             IF ( nr_1d(k) * xrmin > qr_1d(k) * hyrho(k) )  THEN
1424                nr_1d(k) = qr_1d(k) * hyrho(k) / xrmin
1425             ELSEIF ( nr_1d(k) * xrmax < qr_1d(k) * hyrho(k) )  THEN
1426                nr_1d(k) = qr_1d(k) * hyrho(k) / xrmax
1427             ENDIF
1428
1429          ENDIF
1430
1431       ENDDO
1432
1433    END SUBROUTINE adjust_cloud_ij
1434
1435
1436!------------------------------------------------------------------------------!
1437! Description:
1438! ------------
1439!> Autoconversion rate (Seifert and Beheng, 2006). Call for grid point i,j
1440!------------------------------------------------------------------------------!
1441    SUBROUTINE autoconversion_ij( i, j )
1442
1443       USE arrays_3d,                                                          &
1444           ONLY:  diss, dzu
1445
1446       USE cloud_parameters,                                                   &
1447           ONLY:  hyrho
1448
1449       USE control_parameters,                                                 &
1450           ONLY:  rho_surface
1451
1452       USE grid_variables,                                                     &
1453           ONLY:  dx, dy
1454
1455       USE indices,                                                            &
1456           ONLY:  nzb_s_inner, nzt
1457
1458       USE kinds
1459
1460       IMPLICIT NONE
1461
1462       INTEGER(iwp) ::  i                 !<
1463       INTEGER(iwp) ::  j                 !<
1464       INTEGER(iwp) ::  k                 !<
1465
1466       REAL(wp)     ::  alpha_cc          !<                   
1467       REAL(wp)     ::  autocon           !<
1468       REAL(wp)     ::  dissipation       !<
1469       REAL(wp)     ::  k_au              !<
1470       REAL(wp)     ::  l_mix             !<
1471       REAL(wp)     ::  nu_c              !<
1472       REAL(wp)     ::  phi_au            !<
1473       REAL(wp)     ::  r_cc              !<
1474       REAL(wp)     ::  rc                !<
1475       REAL(wp)     ::  re_lambda         !<
1476       REAL(wp)     ::  sigma_cc          !<
1477       REAL(wp)     ::  tau_cloud         !<
1478       REAL(wp)     ::  xc                !<
1479
1480       DO  k = nzb_s_inner(j,i)+1, nzt
1481
1482          IF ( qc_1d(k) > eps_sb )  THEN
1483
1484             k_au = k_cc / ( 20.0_wp * x0 )
1485!
1486!--          Intern time scale of coagulation (Seifert and Beheng, 2006):
1487!--          (1.0_wp - qc(k,j,i) / ( qc(k,j,i) + qr_1d(k) ))
1488             tau_cloud = 1.0_wp - qc_1d(k) / ( qr_1d(k) + qc_1d(k) )
1489!
1490!--          Universal function for autoconversion process
1491!--          (Seifert and Beheng, 2006):
1492             phi_au = 600.0_wp * tau_cloud**0.68_wp * ( 1.0_wp - tau_cloud**0.68_wp )**3
1493!
1494!--          Shape parameter of gamma distribution (Geoffroy et al., 2010):
1495!--          (Use constant nu_c = 1.0_wp instead?)
1496             nu_c = 1.0_wp !MAX( 0.0_wp, 1580.0_wp * hyrho(k) * qc_1d(k) - 0.28_wp )
1497!
1498!--          Mean weight of cloud droplets:
1499             xc = hyrho(k) * qc_1d(k) / nc_1d(k)
1500!
1501!--          Parameterized turbulence effects on autoconversion (Seifert,
1502!--          Nuijens and Stevens, 2010)
1503             IF ( collision_turbulence )  THEN
1504!
1505!--             Weight averaged radius of cloud droplets:
1506                rc = 0.5_wp * ( xc * dpirho_l )**( 1.0_wp / 3.0_wp )
1507
1508                alpha_cc = ( a_1 + a_2 * nu_c ) / ( 1.0_wp + a_3 * nu_c )
1509                r_cc     = ( b_1 + b_2 * nu_c ) / ( 1.0_wp + b_3 * nu_c )
1510                sigma_cc = ( c_1 + c_2 * nu_c ) / ( 1.0_wp + c_3 * nu_c )
1511!
1512!--             Mixing length (neglecting distance to ground and stratification)
1513                l_mix = ( dx * dy * dzu(k) )**( 1.0_wp / 3.0_wp )
1514!
1515!--             Limit dissipation rate according to Seifert, Nuijens and
1516!--             Stevens (2010)
1517                dissipation = MIN( 0.06_wp, diss(k,j,i) )
1518!
1519!--             Compute Taylor-microscale Reynolds number:
1520                re_lambda = 6.0_wp / 11.0_wp *                                 &
1521                            ( l_mix / c_const )**( 2.0_wp / 3.0_wp ) *         &
1522                            SQRT( 15.0_wp / kin_vis_air ) *                    &
1523                            dissipation**( 1.0_wp / 6.0_wp )
1524!
1525!--             The factor of 1.0E4 is needed to convert the dissipation rate
1526!--             from m2 s-3 to cm2 s-3.
1527                k_au = k_au * ( 1.0_wp +                                       &
1528                                dissipation * 1.0E4_wp *                       &
1529                                ( re_lambda * 1.0E-3_wp )**0.25_wp *           &
1530                                ( alpha_cc * EXP( -1.0_wp * ( ( rc - r_cc ) /  &
1531                                                  sigma_cc )**2                &
1532                                                ) + beta_cc                    &
1533                                )                                              &
1534                              )
1535             ENDIF
1536!
1537!--          Autoconversion rate (Seifert and Beheng, 2006):
1538             autocon = k_au * ( nu_c + 2.0_wp ) * ( nu_c + 4.0_wp ) /          &
1539                       ( nu_c + 1.0_wp )**2 * qc_1d(k)**2 * xc**2 *            &
1540                       ( 1.0_wp + phi_au / ( 1.0_wp - tau_cloud )**2 ) *       &
1541                       rho_surface
1542             autocon = MIN( autocon, qc_1d(k) / dt_micro )
1543
1544             qr_1d(k) = qr_1d(k) + autocon * dt_micro
1545             qc_1d(k) = qc_1d(k) - autocon * dt_micro 
1546             nr_1d(k) = nr_1d(k) + autocon / x0 * hyrho(k) * dt_micro
1547
1548          ENDIF
1549
1550       ENDDO
1551
1552    END SUBROUTINE autoconversion_ij
1553
1554!------------------------------------------------------------------------------!
1555! Description:
1556! ------------
1557!> Autoconversion process (Kessler, 1969).
1558!------------------------------------------------------------------------------!
1559    SUBROUTINE autoconversion_kessler_ij( i, j )
1560
1561       USE arrays_3d,                                                          &
1562           ONLY:  dzw, prr
1563
1564       USE cloud_parameters,                                                   &
1565           ONLY:  l_d_cp, pt_d_t
1566
1567       USE indices,                                                            &
1568           ONLY:  nzb_s_inner, nzt
1569
1570       USE kinds
1571
1572
1573       IMPLICIT NONE
1574
1575       INTEGER(iwp) ::  i !<
1576       INTEGER(iwp) ::  j !<
1577       INTEGER(iwp) ::  k !<
1578
1579       REAL(wp)    ::  dqdt_precip !<
1580
1581       DO  k = nzb_s_inner(j,i)+1, nzt
1582
1583          IF ( qc_1d(k) > ql_crit )  THEN
1584             dqdt_precip = prec_time_const * ( qc_1d(k) - ql_crit )
1585          ELSE
1586             dqdt_precip = 0.0_wp
1587          ENDIF
1588
1589          qc_1d(k) = qc_1d(k) - dqdt_precip * dt_micro
1590          q_1d(k)  = q_1d(k)  - dqdt_precip * dt_micro
1591          pt_1d(k) = pt_1d(k) + dqdt_precip * dt_micro * l_d_cp * pt_d_t(k)
1592
1593!
1594!--       Compute the rain rate (stored on surface grid point)
1595          prr(nzb_s_inner(j,i),j,i) = prr(nzb_s_inner(j,i),j,i) +              &
1596                                      dqdt_precip * dzw(k)
1597
1598       ENDDO
1599
1600    END SUBROUTINE autoconversion_kessler_ij
1601
1602!------------------------------------------------------------------------------!
1603! Description:
1604! ------------
1605!> Accretion rate (Seifert and Beheng, 2006). Call for grid point i,j
1606!------------------------------------------------------------------------------!
1607    SUBROUTINE accretion_ij( i, j )
1608
1609       USE arrays_3d,                                                          &
1610           ONLY:  diss
1611
1612       USE cloud_parameters,                                                   &
1613           ONLY:  hyrho
1614
1615       USE control_parameters,                                                 &
1616           ONLY:  rho_surface
1617
1618       USE indices,                                                            &
1619           ONLY:  nzb_s_inner, nzt
1620
1621       USE kinds
1622
1623       IMPLICIT NONE
1624
1625       INTEGER(iwp) ::  i                 !<
1626       INTEGER(iwp) ::  j                 !<
1627       INTEGER(iwp) ::  k                 !<
1628
1629       REAL(wp)     ::  accr              !<
1630       REAL(wp)     ::  k_cr              !<
1631       REAL(wp)     ::  phi_ac            !<
1632       REAL(wp)     ::  tau_cloud         !<
1633
1634       DO  k = nzb_s_inner(j,i)+1, nzt
1635          IF ( ( qc_1d(k) > eps_sb )  .AND.  ( qr_1d(k) > eps_sb ) )  THEN
1636!
1637!--          Intern time scale of coagulation (Seifert and Beheng, 2006):
1638             tau_cloud = 1.0_wp - qc_1d(k) / ( qc_1d(k) + qr_1d(k) ) 
1639!
1640!--          Universal function for accretion process
1641!--          (Seifert and Beheng, 2001):
1642             phi_ac = ( tau_cloud / ( tau_cloud + 5.0E-5_wp ) )**4
1643!
1644!--          Parameterized turbulence effects on autoconversion (Seifert,
1645!--          Nuijens and Stevens, 2010). The factor of 1.0E4 is needed to
1646!--          convert the dissipation rate (diss) from m2 s-3 to cm2 s-3.
1647             IF ( collision_turbulence )  THEN
1648                k_cr = k_cr0 * ( 1.0_wp + 0.05_wp *                            &
1649                                 MIN( 600.0_wp,                                &
1650                                      diss(k,j,i) * 1.0E4_wp )**0.25_wp        &
1651                               )
1652             ELSE
1653                k_cr = k_cr0                       
1654             ENDIF
1655!
1656!--          Accretion rate (Seifert and Beheng, 2006):
1657             accr = k_cr * qc_1d(k) * qr_1d(k) * phi_ac * SQRT( rho_surface * hyrho(k) )
1658             accr = MIN( accr, qc_1d(k) / dt_micro )
1659
1660             qr_1d(k) = qr_1d(k) + accr * dt_micro 
1661             qc_1d(k) = qc_1d(k) - accr * dt_micro
1662
1663          ENDIF
1664
1665       ENDDO
1666
1667    END SUBROUTINE accretion_ij
1668
1669
1670!------------------------------------------------------------------------------!
1671! Description:
1672! ------------
1673!> Collisional breakup rate (Seifert, 2008). Call for grid point i,j
1674!------------------------------------------------------------------------------!
1675    SUBROUTINE selfcollection_breakup_ij( i, j )
1676
1677       USE cloud_parameters,                                                   &
1678           ONLY:  hyrho
1679
1680       USE control_parameters,                                                 &
1681           ONLY:  rho_surface
1682
1683       USE indices,                                                            &
1684           ONLY:  nzb_s_inner, nzt
1685
1686       USE kinds
1687   
1688       IMPLICIT NONE
1689
1690       INTEGER(iwp) ::  i                 !<
1691       INTEGER(iwp) ::  j                 !<
1692       INTEGER(iwp) ::  k                 !<
1693
1694       REAL(wp)     ::  breakup           !<
1695       REAL(wp)     ::  dr                !<
1696       REAL(wp)     ::  phi_br            !<
1697       REAL(wp)     ::  selfcoll          !<
1698
1699       DO  k = nzb_s_inner(j,i)+1, nzt
1700          IF ( qr_1d(k) > eps_sb )  THEN
1701!
1702!--          Selfcollection rate (Seifert and Beheng, 2001):
1703             selfcoll = k_rr * nr_1d(k) * qr_1d(k) * SQRT( hyrho(k) * rho_surface )
1704!
1705!--          Weight averaged diameter of rain drops:
1706             dr = ( hyrho(k) * qr_1d(k) / nr_1d(k) * dpirho_l )**( 1.0_wp / 3.0_wp )
1707!
1708!--          Collisional breakup rate (Seifert, 2008):
1709             IF ( dr >= 0.3E-3_wp )  THEN
1710                phi_br  = k_br * ( dr - 1.1E-3_wp )
1711                breakup = selfcoll * ( phi_br + 1.0_wp )
1712             ELSE
1713                breakup = 0.0_wp
1714             ENDIF
1715
1716             selfcoll = MAX( breakup - selfcoll, -nr_1d(k) / dt_micro )
1717             nr_1d(k) = nr_1d(k) + selfcoll * dt_micro
1718
1719          ENDIF         
1720       ENDDO
1721
1722    END SUBROUTINE selfcollection_breakup_ij
1723
1724
1725!------------------------------------------------------------------------------!
1726! Description:
1727! ------------
1728!> Evaporation of precipitable water. Condensation is neglected for
1729!> precipitable water. Call for grid point i,j
1730!------------------------------------------------------------------------------!
1731    SUBROUTINE evaporation_rain_ij( i, j )
1732
1733       USE arrays_3d,                                                          &
1734           ONLY:  hyp
1735
1736       USE cloud_parameters,                                                   &
1737           ONLY:  hyrho, l_d_cp, l_d_r, l_v, r_v, t_d_pt
1738
1739       USE constants,                                                          &
1740           ONLY:  pi
1741
1742       USE indices,                                                            &
1743           ONLY:  nzb_s_inner, nzt
1744
1745       USE kinds
1746
1747       IMPLICIT NONE
1748
1749       INTEGER(iwp) ::  i                 !<
1750       INTEGER(iwp) ::  j                 !<
1751       INTEGER(iwp) ::  k                 !<
1752
1753       REAL(wp)     ::  alpha             !<
1754       REAL(wp)     ::  dr                !<
1755       REAL(wp)     ::  e_s               !<
1756       REAL(wp)     ::  evap              !<
1757       REAL(wp)     ::  evap_nr           !<
1758       REAL(wp)     ::  f_vent            !<
1759       REAL(wp)     ::  g_evap            !<
1760       REAL(wp)     ::  lambda_r          !<
1761       REAL(wp)     ::  mu_r              !<
1762       REAL(wp)     ::  mu_r_2            !<
1763       REAL(wp)     ::  mu_r_5d2          !<
1764       REAL(wp)     ::  nr_0              !<
1765       REAL(wp)     ::  q_s               !<
1766       REAL(wp)     ::  sat               !<
1767       REAL(wp)     ::  t_l               !<
1768       REAL(wp)     ::  temp              !<
1769       REAL(wp)     ::  xr                !<
1770
1771       DO  k = nzb_s_inner(j,i)+1, nzt
1772          IF ( qr_1d(k) > eps_sb )  THEN
1773!
1774!--          Actual liquid water temperature:
1775             t_l = t_d_pt(k) * pt_1d(k)
1776!
1777!--          Saturation vapor pressure at t_l:
1778             e_s = 610.78_wp * EXP( 17.269_wp * ( t_l - 273.16_wp ) /          &
1779                                    ( t_l - 35.86_wp )                         &
1780                                  )
1781!
1782!--          Computation of saturation humidity:
1783             q_s   = 0.622_wp * e_s / ( hyp(k) - 0.378_wp * e_s )
1784             alpha = 0.622_wp * l_d_r * l_d_cp / ( t_l * t_l )
1785             q_s   = q_s * ( 1.0_wp + alpha * q_1d(k) ) / ( 1.0_wp + alpha * q_s )
1786!
1787!--          Supersaturation:
1788             sat   = ( q_1d(k) - qr_1d(k) - qc_1d(k) ) / q_s - 1.0_wp
1789!
1790!--          Evaporation needs only to be calculated in subsaturated regions
1791             IF ( sat < 0.0_wp )  THEN
1792!
1793!--             Actual temperature:
1794                temp = t_l + l_d_cp * ( qc_1d(k) + qr_1d(k) )
1795       
1796                g_evap = 1.0_wp / ( ( l_v / ( r_v * temp ) - 1.0_wp ) * l_v /  &
1797                                    ( thermal_conductivity_l * temp ) +        &
1798                                    r_v * temp / ( diff_coeff_l * e_s )        &
1799                                  )
1800!
1801!--             Mean weight of rain drops
1802                xr = hyrho(k) * qr_1d(k) / nr_1d(k)
1803!
1804!--             Weight averaged diameter of rain drops:
1805                dr = ( xr * dpirho_l )**( 1.0_wp / 3.0_wp )
1806!
1807!--             Compute ventilation factor and intercept parameter
1808!--             (Seifert and Beheng, 2006; Seifert, 2008):
1809                IF ( ventilation_effect )  THEN
1810!
1811!--                Shape parameter of gamma distribution (Milbrandt and Yau, 2005;
1812!--                Stevens and Seifert, 2008):
1813                   mu_r = 10.0_wp * ( 1.0_wp + TANH( 1.2E3_wp * ( dr - 1.4E-3_wp ) ) )
1814!
1815!--                Slope parameter of gamma distribution (Seifert, 2008):
1816                   lambda_r = ( ( mu_r + 3.0_wp ) * ( mu_r + 2.0_wp ) *        &
1817                                ( mu_r + 1.0_wp )                              &
1818                              )**( 1.0_wp / 3.0_wp ) / dr
1819
1820                   mu_r_2   = mu_r + 2.0_wp
1821                   mu_r_5d2 = mu_r + 2.5_wp 
1822
1823                   f_vent = a_vent * gamm( mu_r_2 ) * lambda_r**( -mu_r_2 ) +  & 
1824                            b_vent * schmidt_p_1d3 *                           &
1825                            SQRT( a_term / kin_vis_air ) * gamm( mu_r_5d2 ) *  &
1826                            lambda_r**( -mu_r_5d2 ) *                          &
1827                            ( 1.0_wp -                                         &
1828                              0.5_wp * ( b_term / a_term ) *                   &
1829                              ( lambda_r / ( c_term + lambda_r )               &
1830                              )**mu_r_5d2 -                                    &
1831                              0.125_wp * ( b_term / a_term )**2 *              &
1832                              ( lambda_r / ( 2.0_wp * c_term + lambda_r )      &
1833                              )**mu_r_5d2 -                                    &
1834                              0.0625_wp * ( b_term / a_term )**3 *             &
1835                              ( lambda_r / ( 3.0_wp * c_term + lambda_r )      &
1836                              )**mu_r_5d2 -                                    &
1837                              0.0390625_wp * ( b_term / a_term )**4 *          & 
1838                              ( lambda_r / ( 4.0_wp * c_term + lambda_r )      &
1839                              )**mu_r_5d2                                      &
1840                            )
1841
1842                   nr_0   = nr_1d(k) * lambda_r**( mu_r + 1.0_wp ) /           &
1843                            gamm( mu_r + 1.0_wp ) 
1844                ELSE
1845                   f_vent = 1.0_wp
1846                   nr_0   = nr_1d(k) * dr
1847                ENDIF
1848!
1849!--             Evaporation rate of rain water content (Seifert and Beheng, 2006):
1850                evap    = 2.0_wp * pi * nr_0 * g_evap * f_vent * sat / hyrho(k)
1851                evap    = MAX( evap, -qr_1d(k) / dt_micro )
1852                evap_nr = MAX( c_evap * evap / xr * hyrho(k),                  &
1853                               -nr_1d(k) / dt_micro )
1854
1855                qr_1d(k) = qr_1d(k) + evap    * dt_micro
1856                nr_1d(k) = nr_1d(k) + evap_nr * dt_micro
1857
1858             ENDIF
1859          ENDIF         
1860
1861       ENDDO
1862
1863    END SUBROUTINE evaporation_rain_ij
1864
1865
1866!------------------------------------------------------------------------------!
1867! Description:
1868! ------------
1869!> Sedimentation of cloud droplets (Ackermann et al., 2009, MWR).
1870!> Call for grid point i,j
1871!------------------------------------------------------------------------------!
1872    SUBROUTINE sedimentation_cloud_ij( i, j )
1873
1874       USE arrays_3d,                                                          &
1875           ONLY:  ddzu, dzu, prr
1876
1877       USE cloud_parameters,                                                   &
1878           ONLY:  hyrho, l_d_cp, pt_d_t
1879
1880       USE control_parameters,                                                 &
1881           ONLY:  call_microphysics_at_all_substeps, intermediate_timestep_count
1882
1883       USE indices,                                                            &
1884           ONLY:  nzb, nzb_s_inner, nzt
1885
1886       USE kinds
1887       
1888       USE statistics,                                                         &
1889           ONLY:  weight_substep
1890
1891       IMPLICIT NONE
1892
1893       INTEGER(iwp) ::  i             !<
1894       INTEGER(iwp) ::  j             !<
1895       INTEGER(iwp) ::  k             !<
1896
1897       REAL(wp), DIMENSION(nzb:nzt+1) :: sed_qc  !<
1898
1899       sed_qc(nzt+1) = 0.0_wp
1900
1901       DO  k = nzt, nzb_s_inner(j,i)+1, -1
1902          IF ( qc_1d(k) > eps_sb )  THEN
1903             sed_qc(k) = sed_qc_const * nc_1d(k)**( -2.0_wp / 3.0_wp ) *       &
1904                         ( qc_1d(k) * hyrho(k) )**( 5.0_wp / 3.0_wp )
1905          ELSE
1906             sed_qc(k) = 0.0_wp
1907          ENDIF
1908
1909          sed_qc(k) = MIN( sed_qc(k), hyrho(k) * dzu(k+1) * q_1d(k) /          &
1910                                      dt_micro + sed_qc(k+1)                   &
1911                         )
1912
1913          q_1d(k)  = q_1d(k)  + ( sed_qc(k+1) - sed_qc(k) ) * ddzu(k+1) /      &
1914                                hyrho(k) * dt_micro
1915          qc_1d(k) = qc_1d(k) + ( sed_qc(k+1) - sed_qc(k) ) * ddzu(k+1) /      & 
1916                                hyrho(k) * dt_micro
1917          pt_1d(k) = pt_1d(k) - ( sed_qc(k+1) - sed_qc(k) ) * ddzu(k+1) /      &
1918                                hyrho(k) * l_d_cp * pt_d_t(k) * dt_micro
1919
1920!
1921!--       Compute the precipitation rate of cloud (fog) droplets
1922          IF ( call_microphysics_at_all_substeps )  THEN
1923             prr(k,j,i) = prr(k,j,i) + sed_qc(k) / hyrho(k) *               &
1924                             weight_substep(intermediate_timestep_count)
1925          ELSE
1926             prr(k,j,i) = prr(k,j,i) + sed_qc(k) / hyrho(k)
1927          ENDIF
1928
1929       ENDDO
1930
1931    END SUBROUTINE sedimentation_cloud_ij
1932
1933
1934!------------------------------------------------------------------------------!
1935! Description:
1936! ------------
1937!> Computation of sedimentation flux. Implementation according to Stevens
1938!> and Seifert (2008). Code is based on UCLA-LES. Call for grid point i,j
1939!------------------------------------------------------------------------------!
1940    SUBROUTINE sedimentation_rain_ij( i, j )
1941
1942       USE arrays_3d,                                                          &
1943           ONLY:  ddzu, dzu, prr
1944
1945       USE cloud_parameters,                                                   &
1946           ONLY:  hyrho, l_d_cp, pt_d_t
1947
1948       USE control_parameters,                                                 &
1949           ONLY:  call_microphysics_at_all_substeps, intermediate_timestep_count
1950
1951       USE indices,                                                            &
1952           ONLY:  nzb, nzb_s_inner, nzt
1953
1954       USE kinds
1955
1956       USE statistics,                                                         &
1957           ONLY:  weight_substep
1958       
1959       IMPLICIT NONE
1960
1961       INTEGER(iwp) ::  i                          !<
1962       INTEGER(iwp) ::  j                          !<
1963       INTEGER(iwp) ::  k                          !<
1964       INTEGER(iwp) ::  k_run                      !<
1965
1966       REAL(wp)     ::  c_run                      !<
1967       REAL(wp)     ::  d_max                      !<
1968       REAL(wp)     ::  d_mean                     !<
1969       REAL(wp)     ::  d_min                      !<
1970       REAL(wp)     ::  dr                         !<
1971       REAL(wp)     ::  flux                       !<
1972       REAL(wp)     ::  lambda_r                   !<
1973       REAL(wp)     ::  mu_r                       !<
1974       REAL(wp)     ::  z_run                      !<
1975
1976       REAL(wp), DIMENSION(nzb:nzt+1) ::  c_nr     !<
1977       REAL(wp), DIMENSION(nzb:nzt+1) ::  c_qr     !<
1978       REAL(wp), DIMENSION(nzb:nzt+1) ::  nr_slope !<
1979       REAL(wp), DIMENSION(nzb:nzt+1) ::  qr_slope !<
1980       REAL(wp), DIMENSION(nzb:nzt+1) ::  sed_nr   !<
1981       REAL(wp), DIMENSION(nzb:nzt+1) ::  sed_qr   !<
1982       REAL(wp), DIMENSION(nzb:nzt+1) ::  w_nr     !<
1983       REAL(wp), DIMENSION(nzb:nzt+1) ::  w_qr     !<
1984
1985!
1986!--    Compute velocities
1987       DO  k = nzb_s_inner(j,i)+1, nzt
1988          IF ( qr_1d(k) > eps_sb )  THEN
1989!
1990!--          Weight averaged diameter of rain drops:
1991             dr = ( hyrho(k) * qr_1d(k) / nr_1d(k) * dpirho_l )**( 1.0_wp / 3.0_wp )
1992!
1993!--          Shape parameter of gamma distribution (Milbrandt and Yau, 2005;
1994!--          Stevens and Seifert, 2008):
1995             mu_r = 10.0_wp * ( 1.0_wp + TANH( 1.2E3_wp * ( dr - 1.4E-3_wp ) ) )
1996!
1997!--          Slope parameter of gamma distribution (Seifert, 2008):
1998             lambda_r = ( ( mu_r + 3.0_wp ) * ( mu_r + 2.0_wp ) *              &
1999                          ( mu_r + 1.0_wp ) )**( 1.0_wp / 3.0_wp ) / dr
2000
2001             w_nr(k) = MAX( 0.1_wp, MIN( 20.0_wp,                              &
2002                                         a_term - b_term * ( 1.0_wp +          &
2003                                            c_term / lambda_r )**( -1.0_wp *   &
2004                                            ( mu_r + 1.0_wp ) )                &
2005                                        )                                      &
2006                          )
2007             w_qr(k) = MAX( 0.1_wp, MIN( 20.0_wp,                              &
2008                                         a_term - b_term * ( 1.0_wp +          &
2009                                            c_term / lambda_r )**( -1.0_wp *   &
2010                                            ( mu_r + 4.0_wp ) )                &
2011                                       )                                       &
2012                          )
2013          ELSE
2014             w_nr(k) = 0.0_wp
2015             w_qr(k) = 0.0_wp
2016          ENDIF
2017       ENDDO
2018!
2019!--    Adjust boundary values
2020       w_nr(nzb_s_inner(j,i)) = w_nr(nzb_s_inner(j,i)+1)
2021       w_qr(nzb_s_inner(j,i)) = w_qr(nzb_s_inner(j,i)+1)
2022       w_nr(nzt+1) = 0.0_wp
2023       w_qr(nzt+1) = 0.0_wp
2024!
2025!--    Compute Courant number
2026       DO  k = nzb_s_inner(j,i)+1, nzt
2027          c_nr(k) = 0.25_wp * ( w_nr(k-1) + 2.0_wp * w_nr(k) + w_nr(k+1) ) *   &
2028                    dt_micro * ddzu(k)
2029          c_qr(k) = 0.25_wp * ( w_qr(k-1) + 2.0_wp * w_qr(k) + w_qr(k+1) ) *   &
2030                    dt_micro * ddzu(k)
2031       ENDDO     
2032!
2033!--    Limit slopes with monotonized centered (MC) limiter (van Leer, 1977):
2034       IF ( limiter_sedimentation )  THEN
2035
2036          DO k = nzb_s_inner(j,i)+1, nzt
2037             d_mean = 0.5_wp * ( qr_1d(k+1) - qr_1d(k-1) )
2038             d_min  = qr_1d(k) - MIN( qr_1d(k+1), qr_1d(k), qr_1d(k-1) )
2039             d_max  = MAX( qr_1d(k+1), qr_1d(k), qr_1d(k-1) ) - qr_1d(k)
2040
2041             qr_slope(k) = SIGN(1.0_wp, d_mean) * MIN ( 2.0_wp * d_min,        &
2042                                                        2.0_wp * d_max,        &
2043                                                        ABS( d_mean ) )
2044
2045             d_mean = 0.5_wp * ( nr_1d(k+1) - nr_1d(k-1) )
2046             d_min  = nr_1d(k) - MIN( nr_1d(k+1), nr_1d(k), nr_1d(k-1) )
2047             d_max  = MAX( nr_1d(k+1), nr_1d(k), nr_1d(k-1) ) - nr_1d(k)
2048
2049             nr_slope(k) = SIGN(1.0_wp, d_mean) * MIN ( 2.0_wp * d_min,        &
2050                                                        2.0_wp * d_max,        &
2051                                                        ABS( d_mean ) )
2052          ENDDO
2053
2054       ELSE
2055
2056          nr_slope = 0.0_wp
2057          qr_slope = 0.0_wp
2058
2059       ENDIF
2060
2061       sed_nr(nzt+1) = 0.0_wp
2062       sed_qr(nzt+1) = 0.0_wp
2063!
2064!--    Compute sedimentation flux
2065       DO  k = nzt, nzb_s_inner(j,i)+1, -1
2066!
2067!--       Sum up all rain drop number densities which contribute to the flux
2068!--       through k-1/2
2069          flux  = 0.0_wp
2070          z_run = 0.0_wp ! height above z(k)
2071          k_run = k
2072          c_run = MIN( 1.0_wp, c_nr(k) )
2073          DO WHILE ( c_run > 0.0_wp  .AND.  k_run <= nzt )
2074             flux  = flux + hyrho(k_run) *                                     &
2075                     ( nr_1d(k_run) + nr_slope(k_run) * ( 1.0_wp - c_run ) *   &
2076                     0.5_wp ) * c_run * dzu(k_run)
2077             z_run = z_run + dzu(k_run)
2078             k_run = k_run + 1
2079             c_run = MIN( 1.0_wp, c_nr(k_run) - z_run * ddzu(k_run) )
2080          ENDDO
2081!
2082!--       It is not allowed to sediment more rain drop number density than
2083!--       available
2084          flux = MIN( flux,                                                    &
2085                      hyrho(k) * dzu(k+1) * nr_1d(k) + sed_nr(k+1) * dt_micro )
2086
2087          sed_nr(k) = flux / dt_micro
2088          nr_1d(k)  = nr_1d(k) + ( sed_nr(k+1) - sed_nr(k) ) * ddzu(k+1) /     &
2089                                    hyrho(k) * dt_micro
2090!
2091!--       Sum up all rain water content which contributes to the flux
2092!--       through k-1/2
2093          flux  = 0.0_wp
2094          z_run = 0.0_wp ! height above z(k)
2095          k_run = k
2096          c_run = MIN( 1.0_wp, c_qr(k) )
2097
2098          DO WHILE ( c_run > 0.0_wp  .AND.  k_run <= nzt )
2099
2100             flux  = flux + hyrho(k_run) *                                     &
2101                     ( qr_1d(k_run) + qr_slope(k_run) * ( 1.0_wp - c_run ) *   &
2102                     0.5_wp ) * c_run * dzu(k_run)
2103             z_run = z_run + dzu(k_run)
2104             k_run = k_run + 1
2105             c_run = MIN( 1.0_wp, c_qr(k_run) - z_run * ddzu(k_run) )
2106
2107          ENDDO
2108!
2109!--       It is not allowed to sediment more rain water content than available
2110          flux = MIN( flux,                                                    &
2111                      hyrho(k) * dzu(k) * qr_1d(k) + sed_qr(k+1) * dt_micro )
2112
2113          sed_qr(k) = flux / dt_micro
2114
2115          qr_1d(k) = qr_1d(k) + ( sed_qr(k+1) - sed_qr(k) ) * ddzu(k+1) /      &
2116                                hyrho(k) * dt_micro
2117          q_1d(k)  = q_1d(k)  + ( sed_qr(k+1) - sed_qr(k) ) * ddzu(k+1) /      &
2118                                hyrho(k) * dt_micro 
2119          pt_1d(k) = pt_1d(k) - ( sed_qr(k+1) - sed_qr(k) ) * ddzu(k+1) /      &
2120                                hyrho(k) * l_d_cp * pt_d_t(k) * dt_micro
2121!
2122!--       Compute the rain rate
2123          IF ( call_microphysics_at_all_substeps )  THEN
2124             prr(k,j,i) = prr(k,j,i) + sed_qr(k) / hyrho(k)                    &
2125                          * weight_substep(intermediate_timestep_count)
2126          ELSE
2127             prr(k,j,i) = prr(k,j,i) + sed_qr(k) / hyrho(k)
2128          ENDIF
2129
2130       ENDDO
2131
2132    END SUBROUTINE sedimentation_rain_ij
2133
2134
2135!------------------------------------------------------------------------------!
2136! Description:
2137! ------------
2138!> This subroutine computes the precipitation amount due to gravitational
2139!> settling of rain and cloud (fog) droplets
2140!------------------------------------------------------------------------------!
2141    SUBROUTINE calc_precipitation_amount_ij( i, j )
2142
2143       USE arrays_3d,                                                          &
2144           ONLY:  precipitation_amount, prr
2145
2146       USE cloud_parameters,                                                   &
2147           ONLY:  hyrho
2148
2149       USE control_parameters,                                                 &
2150           ONLY:  call_microphysics_at_all_substeps, dt_do2d_xy, dt_3d,        &
2151                  intermediate_timestep_count, intermediate_timestep_count_max,&
2152                  precipitation_amount_interval, time_do2d_xy
2153
2154       USE indices,                                                            &
2155           ONLY:  nzb_s_inner
2156
2157       USE kinds
2158
2159       IMPLICIT NONE
2160
2161       INTEGER(iwp) ::  i                          !:
2162       INTEGER(iwp) ::  j                          !:
2163
2164
2165       IF ( ( dt_do2d_xy - time_do2d_xy ) < precipitation_amount_interval .AND.&
2166            ( .NOT. call_microphysics_at_all_substeps .OR.                     &
2167            intermediate_timestep_count == intermediate_timestep_count_max ) ) &
2168       THEN
2169
2170          precipitation_amount(j,i) = precipitation_amount(j,i) +              &
2171                                      prr(nzb_s_inner(j,i)+1,j,i) *            &
2172                                      hyrho(nzb_s_inner(j,i)+1) * dt_3d
2173       ENDIF
2174
2175    END SUBROUTINE calc_precipitation_amount_ij
2176
2177!------------------------------------------------------------------------------!
2178! Description:
2179! ------------
2180!> This function computes the gamma function (Press et al., 1992).
2181!> The gamma function is needed for the calculation of the evaporation
2182!> of rain drops.
2183!------------------------------------------------------------------------------!
2184    FUNCTION gamm( xx ) 
2185
2186       USE kinds
2187
2188       IMPLICIT NONE
2189
2190       INTEGER(iwp) ::  j            !<
2191
2192       REAL(wp)     ::  gamm         !<
2193       REAL(wp)     ::  ser          !<
2194       REAL(wp)     ::  tmp          !<
2195       REAL(wp)     ::  x_gamm       !<
2196       REAL(wp)     ::  xx           !<
2197       REAL(wp)     ::  y_gamm       !<
2198
2199
2200       REAL(wp), PARAMETER  ::  stp = 2.5066282746310005_wp               !<
2201       REAL(wp), PARAMETER  ::  cof(6) = (/ 76.18009172947146_wp,      &
2202                                           -86.50532032941677_wp,      &
2203                                            24.01409824083091_wp,      &
2204                                            -1.231739572450155_wp,     &
2205                                             0.1208650973866179E-2_wp, &
2206                                            -0.5395239384953E-5_wp /)     !<
2207
2208       x_gamm = xx
2209       y_gamm = x_gamm
2210       tmp = x_gamm + 5.5_wp
2211       tmp = ( x_gamm + 0.5_wp ) * LOG( tmp ) - tmp
2212       ser = 1.000000000190015_wp
2213
2214       DO  j = 1, 6 
2215          y_gamm = y_gamm + 1.0_wp 
2216          ser    = ser + cof( j ) / y_gamm 
2217       ENDDO
2218
2219!
2220!--    Until this point the algorithm computes the logarithm of the gamma
2221!--    function. Hence, the exponential function is used. 
2222!       gamm = EXP( tmp + LOG( stp * ser / x_gamm ) )
2223       gamm = EXP( tmp ) * stp * ser / x_gamm 
2224
2225       RETURN
2226 
2227    END FUNCTION gamm 
2228
2229 END MODULE microphysics_mod
Note: See TracBrowser for help on using the repository browser.