source: palm/trunk/SOURCE/lpm_droplet_condensation.f90 @ 2638

Last change on this file since 2638 was 2608, checked in by schwenkel, 6 years ago

Inital revision of diagnostic_quantities_mod allows unified calculation of magnus equation and saturion mixing ratio

  • Property svn:keywords set to Id
File size: 17.3 KB
RevLine 
[1682]1!> @file lpm_droplet_condensation.f90
[2000]2!------------------------------------------------------------------------------!
[1036]3! This file is part of PALM.
4!
[2000]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
[1036]9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2101]17! Copyright 1997-2017 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[849]20! Current revisions:
21! ------------------
[2375]22!
23!
[1891]24! Former revisions:
25! -----------------
26! $Id: lpm_droplet_condensation.f90 2608 2017-11-13 14:04:26Z raasch $
[2608]27! Calculation of magnus equation in external module (diagnostic_quantities_mod).
28!
29! 2375 2017-08-29 14:10:28Z schwenkel
[2375]30! Changed ONLY-dependencies
31!
32! 2312 2017-07-14 20:26:51Z hoffmann
[2312]33! Rosenbrock scheme improved. Gas-kinetic effect added.
[1891]34!
[2001]35! 2000 2016-08-20 18:09:15Z knoop
36! Forced header and separation lines into 80 columns
[2312]37!
[1891]38! 1890 2016-04-22 08:52:11Z hoffmann
[1890]39! Some improvements of the Rosenbrock method. If the Rosenbrock method needs more
40! than 40 iterations to find a sufficient time setp, the model is not aborted.
[2312]41! This might lead to small erros. But they can be assumend as negligible, since
42! the maximum timesetp of the Rosenbrock method after 40 iterations will be
43! smaller than 10^-16 s.
44!
[1872]45! 1871 2016-04-15 11:46:09Z hoffmann
46! Initialization of aerosols added.
47!
[1851]48! 1849 2016-04-08 11:33:18Z hoffmann
[1852]49! Interpolation of supersaturation has been removed because it is not in
50! accordance with the release/depletion of latent heat/water vapor in
[1849]51! interaction_droplets_ptq.
52! Calculation of particle Reynolds number has been corrected.
[1852]53! eps_ros added from modules.
[1849]54!
[1832]55! 1831 2016-04-07 13:15:51Z hoffmann
56! curvature_solution_effects moved to particle_attributes
57!
[1823]58! 1822 2016-04-07 07:49:42Z hoffmann
59! Unused variables removed.
60!
[1683]61! 1682 2015-10-07 23:56:08Z knoop
[2312]62! Code annotations made doxygen readable
63!
[1360]64! 1359 2014-04-11 17:15:14Z hoffmann
[2312]65! New particle structure integrated.
[1360]66! Kind definition added to all floating point numbers.
67!
[1347]68! 1346 2014-03-27 13:18:20Z heinze
[2312]69! Bugfix: REAL constants provided with KIND-attribute especially in call of
[1347]70! intrinsic function like MAX, MIN, SIGN
71!
[1323]72! 1322 2014-03-20 16:38:49Z raasch
73! REAL constants defined as wp-kind
74!
[1321]75! 1320 2014-03-20 08:40:49Z raasch
[1320]76! ONLY-attribute added to USE-statements,
77! kind-parameters added to all INTEGER and REAL declaration statements,
78! kinds are defined in new module kinds,
79! comment fields (!:) to be used for variable explanations added to
80! all variable declaration statements
[1072]81!
[1319]82! 1318 2014-03-17 13:35:16Z raasch
83! module interfaces removed
84!
[1093]85! 1092 2013-02-02 11:24:22Z raasch
86! unused variables removed
87!
[1072]88! 1071 2012-11-29 16:54:55Z franke
[1071]89! Ventilation effect for evaporation of large droplets included
90! Check for unreasonable results included in calculation of Rosenbrock method
91! since physically unlikely results were observed and for the same
92! reason the first internal time step in Rosenbrock method should be < 1.0E02 in
93! case of evaporation
94! Unnecessary calculation of ql_int removed
95! Unnecessary calculations in Rosenbrock method (d2rdt2, drdt_m, dt_ros_last)
96! removed
97! Bugfix: factor in calculation of surface tension changed from 0.00155 to
98! 0.000155
[849]99!
[1037]100! 1036 2012-10-22 13:43:42Z raasch
101! code put under GPL (PALM 3.9)
102!
[850]103! 849 2012-03-15 10:35:09Z raasch
104! initial revision (former part of advec_particles)
[849]105!
[850]106!
[849]107! Description:
108! ------------
[1682]109!> Calculates change in droplet radius by condensation/evaporation, using
110!> either an analytic formula or by numerically integrating the radius growth
111!> equation including curvature and solution effects using Rosenbrocks method
112!> (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
113!> The analytical formula and growth equation follow those given in
114!> Rogers and Yau (A short course in cloud physics, 3rd edition, p. 102/103).
[849]115!------------------------------------------------------------------------------!
[1682]116 SUBROUTINE lpm_droplet_condensation (ip,jp,kp)
[849]117
[2312]118
[1320]119    USE arrays_3d,                                                             &
[2312]120        ONLY:  hyp, pt, q, ql_c, ql_v
[849]121
[1320]122    USE cloud_parameters,                                                      &
[2375]123        ONLY:  l_d_rv, l_v, molecular_weight_of_solute,                        &
124               molecular_weight_of_water, rho_l, rho_s, r_v, vanthoff
[849]125
[1320]126    USE constants,                                                             &
127        ONLY:  pi
[849]128
[1320]129    USE control_parameters,                                                    &
[1822]130        ONLY:  dt_3d, dz, message_string, molecular_viscosity, rho_surface
131
[1320]132    USE cpulog,                                                                &
133        ONLY:  cpu_log, log_point_s
[849]134
[2608]135    USE diagnostic_quantities_mod,                                             &
136        ONLY:  magnus
137
[1320]138    USE grid_variables,                                                        &
[1822]139        ONLY:  dx, dy
[1071]140
[1320]141    USE lpm_collision_kernels_mod,                                             &
142        ONLY:  rclass_lbound, rclass_ubound
[849]143
[1320]144    USE kinds
145
146    USE particle_attributes,                                                   &
[2375]147        ONLY:  curvature_solution_effects, hall_kernel, number_of_particles,   &
148               particles, radius_classes, use_kernel_tables, wang_kernel
[1320]149
150    IMPLICIT NONE
151
[1682]152    INTEGER(iwp) :: ip                         !<
153    INTEGER(iwp) :: jp                         !<
154    INTEGER(iwp) :: kp                         !<
155    INTEGER(iwp) :: n                          !<
[1320]156
[1849]157    REAL(wp) ::  afactor                       !< curvature effects
[1682]158    REAL(wp) ::  arg                           !<
[1849]159    REAL(wp) ::  bfactor                       !< solute effects
[1682]160    REAL(wp) ::  ddenom                        !<
161    REAL(wp) ::  delta_r                       !<
[1849]162    REAL(wp) ::  diameter                      !< diameter of cloud droplets
[2312]163    REAL(wp) ::  diff_coeff                    !< diffusivity for water vapor
[1682]164    REAL(wp) ::  drdt                          !<
165    REAL(wp) ::  dt_ros                        !<
166    REAL(wp) ::  dt_ros_sum                    !<
167    REAL(wp) ::  d2rdtdr                       !<
[1849]168    REAL(wp) ::  e_a                           !< current vapor pressure
169    REAL(wp) ::  e_s                           !< current saturation vapor pressure
[2312]170    REAL(wp) ::  error                         !< local truncation error in Rosenbrock
171    REAL(wp) ::  k1                            !<
172    REAL(wp) ::  k2                            !<
173    REAL(wp) ::  r_err                         !< First order estimate of Rosenbrock radius
174    REAL(wp) ::  r_ros                         !< Rosenbrock radius
175    REAL(wp) ::  r_ros_ini                     !< initial Rosenbrock radius
176    REAL(wp) ::  r0                            !< gas-kinetic lengthscale
177    REAL(wp) ::  sigma                         !< surface tension of water
178    REAL(wp) ::  thermal_conductivity          !< thermal conductivity for water
[1849]179    REAL(wp) ::  t_int                         !< temperature
180    REAL(wp) ::  w_s                           !< terminal velocity of droplets
[2312]181    REAL(wp) ::  re_p                          !< particle Reynolds number
[1849]182!
[2312]183!-- Parameters for Rosenbrock method (see Verwer et al., 1999)
184    REAL(wp), PARAMETER :: prec = 1.0E-3_wp     !< precision of Rosenbrock solution
185    REAL(wp), PARAMETER :: q_increase = 1.5_wp  !< increase factor in timestep
186    REAL(wp), PARAMETER :: q_decrease = 0.9_wp  !< decrease factor in timestep
187    REAL(wp), PARAMETER :: gamma = 0.292893218814_wp !< = 1.0 - 1.0 / SQRT(2.0)
[849]188!
[1849]189!-- Parameters for terminal velocity
190    REAL(wp), PARAMETER ::  a_rog = 9.65_wp      !< parameter for fall velocity
191    REAL(wp), PARAMETER ::  b_rog = 10.43_wp     !< parameter for fall velocity
192    REAL(wp), PARAMETER ::  c_rog = 0.6_wp       !< parameter for fall velocity
193    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp   !< parameter for fall velocity
194    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp  !< parameter for fall velocity
195    REAL(wp), PARAMETER ::  d0_rog = 0.745_wp    !< separation diameter
[849]196
[1849]197    REAL(wp), DIMENSION(number_of_particles) ::  ventilation_effect     !<
198    REAL(wp), DIMENSION(number_of_particles) ::  new_r                  !<
[849]199
[1849]200    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'start' )
[849]201
202!
[2312]203!-- Absolute temperature
[1849]204    t_int = pt(kp,jp,ip) * ( hyp(kp) / 100000.0_wp )**0.286_wp
[849]205!
[2312]206!-- Saturation vapor pressure (Eq. 10 in Bolton, 1980)
[2608]207    e_s = magnus( t_int )
[1849]208!
[2312]209!-- Current vapor pressure
210    e_a = q(kp,jp,ip) * hyp(kp) / ( q(kp,jp,ip) + 0.622_wp )
211!
212!-- Thermal conductivity for water (from Rogers and Yau, Table 7.1)
213    thermal_conductivity = 7.94048E-05_wp * t_int + 0.00227011_wp
214!
215!-- Moldecular diffusivity of water vapor in air (Hall und Pruppacher, 1976)
216    diff_coeff           = 0.211E-4_wp * ( t_int / 273.15_wp )**1.94_wp * &
217                           ( 101325.0_wp / hyp(kp) )
218!
219!-- Lengthscale for gas-kinetic effects (from Mordy, 1959, p. 23):
220    r0 = diff_coeff / 0.036_wp * SQRT( 2.0_wp * pi / ( r_v * t_int ) )
221!
222!-- Calculate effects of heat conductivity and diffusion of water vapor on the
223!-- diffusional growth process (usually known as 1.0 / (F_k + F_d) )
224    ddenom  = 1.0_wp / ( rho_l * r_v * t_int / ( e_s * diff_coeff ) +          &
[1849]225                         ( l_v / ( r_v * t_int ) - 1.0_wp ) * rho_l *          &
[2312]226                         l_v / ( thermal_conductivity * t_int )                &
[1849]227                       )
[1359]228    new_r = 0.0_wp
[1849]229!
230!-- Determine ventilation effect on evaporation of large drops
[1359]231    DO  n = 1, number_of_particles
[1849]232
233       IF ( particles(n)%radius >= 4.0E-5_wp  .AND.  e_a / e_s < 1.0_wp )  THEN
[849]234!
[1849]235!--       Terminal velocity is computed for vertical direction (Rogers et al.,
236!--       1993, J. Appl. Meteorol.)
237          diameter = particles(n)%radius * 2000.0_wp !diameter in mm
238          IF ( diameter <= d0_rog )  THEN
239             w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
240          ELSE
241             w_s = a_rog - b_rog * EXP( -c_rog * diameter )
242          ENDIF
[849]243!
[2312]244!--       Calculate droplet's Reynolds number
[1849]245          re_p = 2.0_wp * particles(n)%radius * w_s / molecular_viscosity
[1071]246!
[1359]247!--       Ventilation coefficient (Rogers and Yau, 1989):
248          IF ( re_p > 2.5_wp )  THEN
[1849]249             ventilation_effect(n) = 0.78_wp + 0.28_wp * SQRT( re_p )
[1071]250          ELSE
[1849]251             ventilation_effect(n) = 1.0_wp + 0.09_wp * re_p
[1071]252          ENDIF
[1849]253       ELSE
[1071]254!
[2312]255!--       For small droplets or in supersaturated environments, the ventilation
[1849]256!--       effect does not play a role
257          ventilation_effect(n) = 1.0_wp
[849]258       ENDIF
[1359]259    ENDDO
[849]260
[2312]261    IF( .NOT. curvature_solution_effects ) then
[849]262!
[2312]263!--    Use analytic model for diffusional growth including gas-kinetic
264!--    effects (Mordy, 1959) but without the impact of aerosols.
[1849]265       DO  n = 1, number_of_particles
[2312]266          arg      = ( particles(n)%radius + r0 )**2 + 2.0_wp * dt_3d * ddenom * &
267                                                       ventilation_effect(n) *   &
268                                                       ( e_a / e_s - 1.0_wp )
269          arg      = MAX( arg, ( 0.01E-6 + r0 )**2 )
270          new_r(n) = SQRT( arg ) - r0
[1849]271       ENDDO
[1359]272
[2312]273    ELSE
[1849]274!
[2312]275!--    Integrate the diffusional growth including gas-kinetic (Mordy, 1959),
276!--    as well as curvature and solute effects (e.g., Köhler, 1936).
[849]277!
[2312]278!--    Curvature effect (afactor) with surface tension (sigma) by Straka (2009)
279       sigma = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
[1071]280!
[2312]281!--    Solute effect (afactor)
282       afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
[849]283
[2312]284       DO  n = 1, number_of_particles
[849]285!
[2312]286!--       Solute effect (bfactor)
287          bfactor = vanthoff * rho_s * particles(n)%aux1**3 *                    &
288                    molecular_weight_of_water / ( rho_l * molecular_weight_of_solute )
[1071]289
[2312]290          dt_ros     = particles(n)%aux2  ! use previously stored Rosenbrock timestep
291          dt_ros_sum = 0.0_wp
[1871]292
[2312]293          r_ros     = particles(n)%radius  ! initialize Rosenbrock particle radius
294          r_ros_ini = r_ros
[849]295!
[2312]296!--       Integrate growth equation using a 2nd-order Rosenbrock method
297!--       (see Verwer et al., 1999, Eq. (3.2)). The Rosenbrock method adjusts
298!--       its with internal timestep to minimize the local truncation error.
299          DO WHILE ( dt_ros_sum < dt_3d )
[1071]300
[2312]301             dt_ros = MIN( dt_ros, dt_3d - dt_ros_sum )
[1871]302
[2312]303             DO
[849]304
[2312]305                drdt = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0 -    &
[1849]306                                                          afactor / r_ros +    &
307                                                          bfactor / r_ros**3   &
[2312]308                                                        ) / ( r_ros + r0 )
[1849]309
[2312]310                d2rdtdr = -ddenom * ventilation_effect(n) * (                  &
311                                                (e_a / e_s - 1.0) * r_ros**4 - &
312                                                afactor * r0 * r_ros**2 -      &
313                                                2.0 * afactor * r_ros**3 +     &
314                                                3.0 * bfactor * r0 +           &
315                                                4.0 * bfactor * r_ros          &
316                                                            )                  &
317                          / ( r_ros**4 * ( r_ros + r0 )**2 )
[849]318
[2312]319                k1    = drdt / ( 1.0 - gamma * dt_ros * d2rdtdr )
[849]320
[2312]321                r_ros = MAX(r_ros_ini + k1 * dt_ros, particles(n)%aux1)
322                r_err = r_ros
[849]323
[2312]324                drdt  = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0 -   &
325                                                           afactor / r_ros +   &
326                                                           bfactor / r_ros**3  &
327                                                         ) / ( r_ros + r0 )
[849]328
[2312]329                k2 = ( drdt - dt_ros * 2.0 * gamma * d2rdtdr * k1 ) / &
330                     ( 1.0 - dt_ros * gamma * d2rdtdr )
[849]331
[2312]332                r_ros = MAX(r_ros_ini + dt_ros * ( 1.5 * k1 + 0.5 * k2), particles(n)%aux1)
333   !
334   !--          Check error of the solution, and reduce dt_ros if necessary.
335                error = ABS(r_err - r_ros) / r_ros
336                IF ( error .GT. prec )  THEN
337                   dt_ros = SQRT( q_decrease * prec / error ) * dt_ros
338                   r_ros  = r_ros_ini
[849]339                ELSE
[2312]340                   dt_ros_sum = dt_ros_sum + dt_ros
341                   dt_ros     = q_increase * dt_ros
342                   r_ros_ini  = r_ros
343                   EXIT
[849]344                ENDIF
345
[2312]346             END DO
[849]347
[2312]348          END DO !Rosenbrock loop
[849]349!
[2312]350!--       Store new particle radius
351          new_r(n) = r_ros
[849]352!
[2312]353!--       Store internal time step value for next PALM step
354          particles(n)%aux2 = dt_ros
[849]355
[2312]356       ENDDO !Particle loop
[849]357
[2312]358    ENDIF
[849]359
[2312]360    DO  n = 1, number_of_particles
[849]361!
[2312]362!--    Sum up the change in liquid water for the respective grid
363!--    box for the computation of the release/depletion of water vapor
364!--    and heat.
[1890]365       ql_c(kp,jp,ip) = ql_c(kp,jp,ip) + particles(n)%weight_factor *          &
[1359]366                                   rho_l * 1.33333333_wp * pi *                &
367                                   ( new_r(n)**3 - particles(n)%radius**3 ) /  &
[849]368                                   ( rho_surface * dx * dy * dz )
[2312]369!
370!--    Check if the increase in liqid water is not too big. If this is the case,
371!--    the model timestep might be too long.
[1890]372       IF ( ql_c(kp,jp,ip) > 100.0_wp )  THEN
373          WRITE( message_string, * ) 'k=',kp,' j=',jp,' i=',ip,      &
374                       ' ql_c=',ql_c(kp,jp,ip), ' &part(',n,')%wf=', &
[849]375                       particles(n)%weight_factor,' delta_r=',delta_r
376          CALL message( 'lpm_droplet_condensation', 'PA0143', 2, 2, -1, 6, 1 )
377       ENDIF
378!
[2312]379!--    Check if the change in the droplet radius is not too big. If this is the
380!--    case, the model timestep might be too long.
381       delta_r = new_r(n) - particles(n)%radius
382       IF ( delta_r < 0.0_wp  .AND. new_r(n) < 0.0_wp )  THEN
[1890]383          WRITE( message_string, * ) '#1 k=',kp,' j=',jp,' i=',ip,    &
[1849]384                       ' e_s=',e_s, ' e_a=',e_a,' t_int=',t_int,      &
[1890]385                       ' &delta_r=',delta_r,                          &
[849]386                       ' particle_radius=',particles(n)%radius
387          CALL message( 'lpm_droplet_condensation', 'PA0144', 2, 2, -1, 6, 1 )
388       ENDIF
389!
390!--    Sum up the total volume of liquid water (needed below for
391!--    re-calculating the weighting factors)
[1890]392       ql_v(kp,jp,ip) = ql_v(kp,jp,ip) + particles(n)%weight_factor * new_r(n)**3
[849]393!
394!--    Determine radius class of the particle needed for collision
[2312]395       IF ( use_kernel_tables )  THEN
[1359]396          particles(n)%class = ( LOG( new_r(n) ) - rclass_lbound ) /           &
397                               ( rclass_ubound - rclass_lbound ) *             &
[849]398                               radius_classes
399          particles(n)%class = MIN( particles(n)%class, radius_classes )
400          particles(n)%class = MAX( particles(n)%class, 1 )
401       ENDIF
[2312]402 !
403 !--   Store new radius to particle features
404       particles(n)%radius = new_r(n)
[849]405
406    ENDDO
407
408    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'stop' )
409
410
411 END SUBROUTINE lpm_droplet_condensation
Note: See TracBrowser for help on using the repository browser.