source: palm/trunk/SOURCE/microphysics.f90 @ 1845

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

nzb_2d replaced by nzb_..._inner, Kessler precipitation stored on surface grid point

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