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

Last change on this file since 1874 was 1872, checked in by hoffmann, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 22.6 KB
RevLine 
[1682]1!> @file lpm_droplet_condensation.f90
[1036]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!
[1818]16! Copyright 1997-2016 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[849]19! Current revisions:
20! ------------------
[1872]21!
[1851]22!
23! Former revisions:
24! -----------------
25! $Id: lpm_droplet_condensation.f90 1872 2016-04-15 11:48:48Z maronga $
[1852]26!
[1872]27! 1871 2016-04-15 11:46:09Z hoffmann
28! Initialization of aerosols added.
29!
[1851]30! 1849 2016-04-08 11:33:18Z hoffmann
[1852]31! Interpolation of supersaturation has been removed because it is not in
32! accordance with the release/depletion of latent heat/water vapor in
[1849]33! interaction_droplets_ptq.
34! Calculation of particle Reynolds number has been corrected.
[1852]35! eps_ros added from modules.
[1849]36!
[1832]37! 1831 2016-04-07 13:15:51Z hoffmann
38! curvature_solution_effects moved to particle_attributes
39!
[1823]40! 1822 2016-04-07 07:49:42Z hoffmann
41! Unused variables removed.
42!
[1683]43! 1682 2015-10-07 23:56:08Z knoop
44! Code annotations made doxygen readable
45!
[1360]46! 1359 2014-04-11 17:15:14Z hoffmann
47! New particle structure integrated.
48! Kind definition added to all floating point numbers.
49!
[1347]50! 1346 2014-03-27 13:18:20Z heinze
51! Bugfix: REAL constants provided with KIND-attribute especially in call of
52! intrinsic function like MAX, MIN, SIGN
53!
[1323]54! 1322 2014-03-20 16:38:49Z raasch
55! REAL constants defined as wp-kind
56!
[1321]57! 1320 2014-03-20 08:40:49Z raasch
[1320]58! ONLY-attribute added to USE-statements,
59! kind-parameters added to all INTEGER and REAL declaration statements,
60! kinds are defined in new module kinds,
61! comment fields (!:) to be used for variable explanations added to
62! all variable declaration statements
[1072]63!
[1319]64! 1318 2014-03-17 13:35:16Z raasch
65! module interfaces removed
66!
[1093]67! 1092 2013-02-02 11:24:22Z raasch
68! unused variables removed
69!
[1072]70! 1071 2012-11-29 16:54:55Z franke
[1071]71! Ventilation effect for evaporation of large droplets included
72! Check for unreasonable results included in calculation of Rosenbrock method
73! since physically unlikely results were observed and for the same
74! reason the first internal time step in Rosenbrock method should be < 1.0E02 in
75! case of evaporation
76! Unnecessary calculation of ql_int removed
77! Unnecessary calculations in Rosenbrock method (d2rdt2, drdt_m, dt_ros_last)
78! removed
79! Bugfix: factor in calculation of surface tension changed from 0.00155 to
80! 0.000155
[849]81!
[1037]82! 1036 2012-10-22 13:43:42Z raasch
83! code put under GPL (PALM 3.9)
84!
[850]85! 849 2012-03-15 10:35:09Z raasch
86! initial revision (former part of advec_particles)
[849]87!
[850]88!
[849]89! Description:
90! ------------
[1682]91!> Calculates change in droplet radius by condensation/evaporation, using
92!> either an analytic formula or by numerically integrating the radius growth
93!> equation including curvature and solution effects using Rosenbrocks method
94!> (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
95!> The analytical formula and growth equation follow those given in
96!> Rogers and Yau (A short course in cloud physics, 3rd edition, p. 102/103).
[849]97!------------------------------------------------------------------------------!
[1682]98 SUBROUTINE lpm_droplet_condensation (ip,jp,kp)
99 
[849]100
[1320]101    USE arrays_3d,                                                             &
[1849]102        ONLY:  hyp, pt, q,  ql_c, ql_v
[849]103
[1320]104    USE cloud_parameters,                                                      &
[1849]105        ONLY:  l_d_rv, l_v, rho_l, r_v
[849]106
[1320]107    USE constants,                                                             &
108        ONLY:  pi
[849]109
[1320]110    USE control_parameters,                                                    &
[1822]111        ONLY:  dt_3d, dz, message_string, molecular_viscosity, rho_surface
112
[1320]113    USE cpulog,                                                                &
114        ONLY:  cpu_log, log_point_s
[849]115
[1320]116    USE grid_variables,                                                        &
[1822]117        ONLY:  dx, dy
[1071]118
[1320]119    USE lpm_collision_kernels_mod,                                             &
120        ONLY:  rclass_lbound, rclass_ubound
[849]121
[1320]122    USE kinds
123
124    USE particle_attributes,                                                   &
[1871]125        ONLY:  curvature_solution_effects, hall_kernel,                        &
[1849]126               molecular_weight_of_solute, molecular_weight_of_water,          &
[1871]127               number_of_particles, particles, radius_classes, rho_s,          &
[1849]128               use_kernel_tables, vanthoff, wang_kernel
[1320]129
130
131    IMPLICIT NONE
132
[1682]133    INTEGER(iwp) :: i                          !<
134    INTEGER(iwp) :: ip                         !<
135    INTEGER(iwp) :: internal_timestep_count    !<
136    INTEGER(iwp) :: j                          !<
137    INTEGER(iwp) :: jp                         !<
138    INTEGER(iwp) :: jtry                       !<
139    INTEGER(iwp) :: k                          !<
140    INTEGER(iwp) :: kp                         !<
141    INTEGER(iwp) :: n                          !<
142    INTEGER(iwp) :: ros_count                  !<
[1320]143 
[1871]144    INTEGER(iwp), PARAMETER ::  maxtry = 100   !<
[1320]145
[1849]146    LOGICAL ::  repeat                         !<
[1320]147
[1682]148    REAL(wp) ::  aa                            !<
[1849]149    REAL(wp) ::  afactor                       !< curvature effects
[1682]150    REAL(wp) ::  arg                           !<
[1849]151    REAL(wp) ::  bfactor                       !< solute effects
[1682]152    REAL(wp) ::  ddenom                        !<
153    REAL(wp) ::  delta_r                       !<
[1849]154    REAL(wp) ::  diameter                      !< diameter of cloud droplets
155    REAL(wp) ::  diff_coeff_v                  !< diffusivity for water vapor
[1682]156    REAL(wp) ::  drdt                          !<
157    REAL(wp) ::  drdt_ini                      !<
158    REAL(wp) ::  dt_ros                        !<
159    REAL(wp) ::  dt_ros_next                   !<
160    REAL(wp) ::  dt_ros_sum                    !<
161    REAL(wp) ::  dt_ros_sum_ini                !<
162    REAL(wp) ::  d2rdtdr                       !<
163    REAL(wp) ::  errmax                        !<
[1849]164    REAL(wp) ::  e_a                           !< current vapor pressure
165    REAL(wp) ::  e_s                           !< current saturation vapor pressure
[1682]166    REAL(wp) ::  err_ros                       !<
167    REAL(wp) ::  g1                            !<
168    REAL(wp) ::  g2                            !<
169    REAL(wp) ::  g3                            !<
170    REAL(wp) ::  g4                            !<
171    REAL(wp) ::  r_ros                         !<
172    REAL(wp) ::  r_ros_ini                     !<
173    REAL(wp) ::  sigma                         !<
[1849]174    REAL(wp) ::  thermal_conductivity_v        !< thermal conductivity for water
175    REAL(wp) ::  t_int                         !< temperature
176    REAL(wp) ::  w_s                           !< terminal velocity of droplets
[1682]177    REAL(wp) ::  re_p                          !<
[1849]178
179!
[849]180!-- Parameters for Rosenbrock method
[1682]181    REAL(wp), PARAMETER ::  a21 = 2.0_wp               !<
182    REAL(wp), PARAMETER ::  a31 = 48.0_wp / 25.0_wp    !<
183    REAL(wp), PARAMETER ::  a32 = 6.0_wp / 25.0_wp     !<
184    REAL(wp), PARAMETER ::  b1 = 19.0_wp / 9.0_wp      !<
185    REAL(wp), PARAMETER ::  b2 = 0.5_wp                !<
186    REAL(wp), PARAMETER ::  b3 = 25.0_wp / 108.0_wp    !<
187    REAL(wp), PARAMETER ::  b4 = 125.0_wp / 108.0_wp   !<
188    REAL(wp), PARAMETER ::  c21 = -8.0_wp              !<
189    REAL(wp), PARAMETER ::  c31 = 372.0_wp / 25.0_wp   !<
190    REAL(wp), PARAMETER ::  c32 = 12.0_wp / 5.0_wp     !<
191    REAL(wp), PARAMETER ::  c41 = -112.0_wp / 125.0_wp !<
192    REAL(wp), PARAMETER ::  c42 = -54.0_wp / 125.0_wp  !<
193    REAL(wp), PARAMETER ::  c43 = -2.0_wp / 5.0_wp     !<
194    REAL(wp), PARAMETER ::  errcon = 0.1296_wp         !<
195    REAL(wp), PARAMETER ::  e1 = 17.0_wp / 54.0_wp     !<
196    REAL(wp), PARAMETER ::  e2 = 7.0_wp / 36.0_wp      !<
197    REAL(wp), PARAMETER ::  e3 = 0.0_wp                !<
198    REAL(wp), PARAMETER ::  e4 = 125.0_wp / 108.0_wp   !<
[1871]199    REAL(wp), PARAMETER ::  eps_ros = 1.0E-3_wp        !< accuracy of Rosenbrock method
[1682]200    REAL(wp), PARAMETER ::  gam = 0.5_wp               !<
201    REAL(wp), PARAMETER ::  grow = 1.5_wp              !<
202    REAL(wp), PARAMETER ::  pgrow = -0.25_wp           !<
203    REAL(wp), PARAMETER ::  pshrnk = -1.0_wp /3.0_wp   !<
204    REAL(wp), PARAMETER ::  shrnk = 0.5_wp             !<
[849]205
206!
[1849]207!-- Parameters for terminal velocity
208    REAL(wp), PARAMETER ::  a_rog = 9.65_wp      !< parameter for fall velocity
209    REAL(wp), PARAMETER ::  b_rog = 10.43_wp     !< parameter for fall velocity
210    REAL(wp), PARAMETER ::  c_rog = 0.6_wp       !< parameter for fall velocity
211    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp   !< parameter for fall velocity
212    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp  !< parameter for fall velocity
213    REAL(wp), PARAMETER ::  d0_rog = 0.745_wp    !< separation diameter
[849]214
[1849]215    REAL(wp), DIMENSION(number_of_particles) ::  ventilation_effect     !<
216    REAL(wp), DIMENSION(number_of_particles) ::  new_r                  !<
[849]217
218
219
[1849]220    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'start' )
[849]221
222!
[1849]223!-- Calculate temperature, saturation vapor pressure and current vapor pressure
224    t_int = pt(kp,jp,ip) * ( hyp(kp) / 100000.0_wp )**0.286_wp
225    e_s   = 611.0_wp * EXP( l_d_rv * ( 3.6609E-3_wp - 1.0_wp / t_int ) )
226    e_a   = q(kp,jp,ip) * hyp(kp) / ( 0.378_wp * q(kp,jp,ip) + 0.622_wp )
[849]227!
[1849]228!-- Thermal conductivity for water (from Rogers and Yau, Table 7.1),
229!-- diffusivity for water vapor (after Hall und Pruppacher, 1976)
230    thermal_conductivity_v = 7.94048E-05_wp * t_int + 0.00227011_wp
231    diff_coeff_v           = 0.211E-4_wp * ( t_int / 273.15_wp )**1.94_wp * &
232                             ( 101325.0_wp / hyp(kp) )
233!
234!-- Calculate effects of heat conductivity and diffusion of water vapor on the
235!-- condensation/evaporation process (typically known as 1.0 / (F_k + F_d) )
236    ddenom  = 1.0_wp / ( rho_l * r_v * t_int / ( e_s * diff_coeff_v ) +        &
237                         ( l_v / ( r_v * t_int ) - 1.0_wp ) * rho_l *          &
238                         l_v / ( thermal_conductivity_v * t_int )              &
239                       )
[849]240
[1359]241    new_r = 0.0_wp
242
[1849]243!
244!-- Determine ventilation effect on evaporation of large drops
[1359]245    DO  n = 1, number_of_particles
[1849]246
247       IF ( particles(n)%radius >= 4.0E-5_wp  .AND.  e_a / e_s < 1.0_wp )  THEN
[849]248!
[1849]249!--       Terminal velocity is computed for vertical direction (Rogers et al.,
250!--       1993, J. Appl. Meteorol.)
251          diameter = particles(n)%radius * 2000.0_wp !diameter in mm
252          IF ( diameter <= d0_rog )  THEN
253             w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
254          ELSE
255             w_s = a_rog - b_rog * EXP( -c_rog * diameter )
256          ENDIF
[849]257!
[1849]258!--       First calculate droplet's Reynolds number
259          re_p = 2.0_wp * particles(n)%radius * w_s / molecular_viscosity
[1071]260!
[1359]261!--       Ventilation coefficient (Rogers and Yau, 1989):
262          IF ( re_p > 2.5_wp )  THEN
[1849]263             ventilation_effect(n) = 0.78_wp + 0.28_wp * SQRT( re_p )
[1071]264          ELSE
[1849]265             ventilation_effect(n) = 1.0_wp + 0.09_wp * re_p
[1071]266          ENDIF
[1849]267       ELSE
[1071]268!
[1849]269!--       For small droplets or in supersaturated environments, the ventilation
270!--       effect does not play a role
271          ventilation_effect(n) = 1.0_wp
[849]272       ENDIF
[1359]273    ENDDO
[849]274
275!
[1849]276!-- Use analytic model for condensational growth
277    IF( .NOT. curvature_solution_effects ) then
278       DO  n = 1, number_of_particles
279          arg = particles(n)%radius**2 + 2.0_wp * dt_3d * ddenom *             &
280                                         ventilation_effect(n) *               &
281                                         ( e_a / e_s - 1.0_wp )
[1359]282          arg = MAX( arg, 1.0E-16_wp )
283          new_r(n) = SQRT( arg )
[1849]284       ENDDO
285    ENDIF
[1359]286
[1849]287!
288!-- If selected, use numerical solution of the condensational growth
289!-- equation (e.g., for studying the activation of aerosols).
290!-- Curvature and solutions effects are included in growth equation.
291!-- Change in Radius is calculated with a 4th-order Rosenbrock method
292!-- for stiff o.d.e's with monitoring local truncation error to adjust
293!-- stepsize (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
[1359]294    DO  n = 1, number_of_particles
[1849]295       IF ( curvature_solution_effects )  THEN
[1071]296
297          ros_count = 0
298          repeat = .TRUE.
[849]299!
[1071]300!--       Carry out the Rosenbrock algorithm. In case of unreasonable results
301!--       the switch "repeat" will be set true and the algorithm will be carried
302!--       out again with the internal time step set to its initial (small) value.
[1359]303!--       Unreasonable results may occur if the external conditions, especially
304!--       the supersaturation, has significantly changed compared to the last
305!--       PALM timestep.
[1071]306          DO WHILE ( repeat )
[849]307
[1071]308             repeat = .FALSE.
309!
[1849]310!--          Curvature effect (afactor) with surface tension parameterization
311!--          by Straka (2009)
312             sigma = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
313             afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
314!
[1871]315!--          Solute effect (bfactor)
316             bfactor = vanthoff * rho_s * particles(n)%rvar2**3 *              &
317                       molecular_weight_of_water /                             &
318                       ( rho_l * molecular_weight_of_solute )
[849]319
[1071]320             r_ros = particles(n)%radius
[1359]321             dt_ros_sum  = 0.0_wp      ! internal integrated time (s)
[1071]322             internal_timestep_count = 0
[849]323!
[1071]324!--          Take internal time step values from the end of last PALM time step
325             dt_ros_next = particles(n)%rvar1
326
[849]327!
[1871]328!--          Internal time step should not be > 1.0E-2 and < 1.0E-6
329             IF ( dt_ros_next > 1.0E-2_wp )  THEN
[1359]330                dt_ros_next = 1.0E-3_wp
[1871]331             ELSEIF ( dt_ros_next < 1.0E-6_wp )  THEN
332                dt_ros_next = 1.0E-6_wp
[1071]333             ENDIF
[1871]334
[849]335!
[1071]336!--          If calculation of Rosenbrock method is repeated due to unreasonalble
337!--          results during previous try the initial internal time step has to be
338!--          reduced
339             IF ( ros_count > 1 )  THEN
[1871]340                dt_ros_next = dt_ros_next * 0.1_wp
[1071]341             ELSEIF ( ros_count > 5 )  THEN
[849]342!
[1071]343!--             Prevent creation of infinite loop
344                message_string = 'ros_count > 5 in Rosenbrock method'
345                CALL message( 'lpm_droplet_condensation', 'PA0018', 2, 2, &
346                               0, 6, 0 )
347             ENDIF
348
[849]349!
[1071]350!--          Internal time step must not be larger than PALM time step
351             dt_ros = MIN( dt_ros_next, dt_3d )
[1871]352
[1071]353!
354!--          Integrate growth equation in time unless PALM time step is reached
355             DO WHILE ( dt_ros_sum < dt_3d )
[849]356
[1071]357                internal_timestep_count = internal_timestep_count + 1
[849]358
359!
[1071]360!--             Derivative at starting value
[1849]361                drdt = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0_wp - &
362                                                          afactor / r_ros +    &
363                                                          bfactor / r_ros**3   &
364                                                        ) / r_ros
365
[1071]366                drdt_ini       = drdt
367                dt_ros_sum_ini = dt_ros_sum
[1871]368                r_ros_ini      = MAX( r_ros, particles(n)%rvar2 )
[849]369
370!
[1071]371!--             Calculate radial derivative of dr/dt
[1849]372                d2rdtdr = ddenom * ventilation_effect(n) *                     &
373                                       ( ( 1.0_wp - e_a / e_s ) / r_ros**2 +   &
374                                         2.0_wp * afactor / r_ros**3 -         &
375                                         4.0_wp * bfactor / r_ros**5           &
376                                       )
[849]377!
[1071]378!--             Adjust stepsize unless required accuracy is reached
379                DO  jtry = 1, maxtry+1
[849]380
[1071]381                   IF ( jtry == maxtry+1 )  THEN
[1871]382                      message_string = 'maxtry > 100 in Rosenbrock method'
[1359]383                      CALL message( 'lpm_droplet_condensation', 'PA0347', 2,   &
384                                    2, 0, 6, 0 )
[1071]385                   ENDIF
[849]386
[1359]387                   aa    = 1.0_wp / ( gam * dt_ros ) - d2rdtdr
[1071]388                   g1    = drdt_ini / aa
[1871]389                   r_ros = MAX( r_ros_ini + a21 * g1, particles(n)%rvar2 )
[1849]390                   drdt  = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0_wp - &
391                                                              afactor / r_ros +    &
392                                                              bfactor / r_ros**3   &
393                                                            ) / r_ros
[849]394
[1071]395                   g2    = ( drdt + c21 * g1 / dt_ros )&
396                           / aa
[1871]397                   r_ros = MAX( r_ros_ini + a31 * g1 + a32 * g2, particles(n)%rvar2 )
[1849]398                   drdt  = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0_wp - &
399                                                              afactor / r_ros +    &
400                                                              bfactor / r_ros**3   &
401                                                            ) / r_ros
[849]402
[1071]403                   g3    = ( drdt +  &
404                             ( c31 * g1 + c32 * g2 ) / dt_ros ) / aa
405                   g4    = ( drdt +  &
406                             ( c41 * g1 + c42 * g2 + c43 * g3 ) / dt_ros ) / aa
[1871]407                   r_ros = MAX( r_ros_ini + b1 * g1 + b2 * g2 + b3 * g3 +      &
408                                b4 * g4, particles(n)%rvar2 )
[849]409
[1071]410                   dt_ros_sum = dt_ros_sum_ini + dt_ros
[849]411
[1071]412                   IF ( dt_ros_sum == dt_ros_sum_ini )  THEN
413                      message_string = 'zero stepsize in Rosenbrock method'
[1359]414                      CALL message( 'lpm_droplet_condensation', 'PA0348', 2,   &
415                                    2, 0, 6, 0 )
[1071]416                   ENDIF
[849]417!
[1071]418!--                Calculate error
[1359]419                   err_ros = e1 * g1 + e2 * g2 + e3 * g3 + e4 * g4
420                   errmax  = 0.0_wp
[1071]421                   errmax  = MAX( errmax, ABS( err_ros / r_ros_ini ) ) / eps_ros
[849]422!
[1071]423!--                Leave loop if accuracy is sufficient, otherwise try again
424!--                with a reduced stepsize
[1359]425                   IF ( errmax <= 1.0_wp )  THEN
[1071]426                      EXIT
427                   ELSE
[1359]428                      dt_ros = SIGN( MAX( ABS( 0.9_wp * dt_ros *               &
429                                               errmax**pshrnk ),               &
430                                               shrnk * ABS( dt_ros ) ), dt_ros )
[1071]431                   ENDIF
432
433                ENDDO  ! loop for stepsize adjustment
434
435!
436!--             Calculate next internal time step
437                IF ( errmax > errcon )  THEN
[1359]438                   dt_ros_next = 0.9_wp * dt_ros * errmax**pgrow
[849]439                ELSE
[1071]440                   dt_ros_next = grow * dt_ros
[849]441                ENDIF
442
[1071]443!
444!--             Estimated time step is reduced if the PALM time step is exceeded
445                IF ( ( dt_ros_next + dt_ros_sum ) >= dt_3d )  THEN
446                   dt_ros = dt_3d - dt_ros_sum
447                ELSE
448                   dt_ros = dt_ros_next
449                ENDIF
[849]450
[1071]451             ENDDO
[849]452!
[1071]453!--          Store internal time step value for next PALM step
454             particles(n)%rvar1 = dt_ros_next
[849]455
456!
[1071]457!--          Radius should not fall below 1E-8 because Rosenbrock method may
458!--          lead to errors otherwise
[1871]459             new_r(n) = MAX( r_ros, particles(n)%rvar2 )
[1071]460!
461!--          Check if calculated droplet radius change is reasonable since in
462!--          case of droplet evaporation the Rosenbrock method may lead to
463!--          secondary solutions which are physically unlikely.
464!--          Due to the solution effect the droplets may grow for relative
[1359]465!--          humidities below 100%, but change of radius should not be too
466!--          large. In case of unreasonable droplet growth the Rosenbrock
467!--          method is recalculated using a smaller initial time step.
[1071]468!--          Limiting values are tested for droplets down to 1.0E-7
[1359]469             IF ( new_r(n) - particles(n)%radius >= 3.0E-7_wp  .AND.  &
[1849]470                  e_a / e_s < 0.97_wp )  THEN
[1071]471                ros_count = ros_count + 1
472                repeat = .TRUE.
[849]473             ENDIF
474
[1071]475          ENDDO    ! Rosenbrock method
[849]476
477       ENDIF
478
[1359]479       delta_r = new_r(n) - particles(n)%radius
[849]480
481!
482!--    Sum up the change in volume of liquid water for the respective grid
483!--    volume (this is needed later in lpm_calc_liquid_water_content for
484!--    calculating the release of latent heat)
[1359]485       i = ip
486       j = jp
487       k = kp
[849]488           ! only exact if equidistant
489
[1359]490       ql_c(k,j,i) = ql_c(k,j,i) + particles(n)%weight_factor *                &
491                                   rho_l * 1.33333333_wp * pi *                &
492                                   ( new_r(n)**3 - particles(n)%radius**3 ) /  &
[849]493                                   ( rho_surface * dx * dy * dz )
[1359]494       IF ( ql_c(k,j,i) > 100.0_wp )  THEN
[849]495          WRITE( message_string, * ) 'k=',k,' j=',j,' i=',i,      &
496                       ' ql_c=',ql_c(k,j,i), ' &part(',n,')%wf=', &
497                       particles(n)%weight_factor,' delta_r=',delta_r
498          CALL message( 'lpm_droplet_condensation', 'PA0143', 2, 2, -1, 6, 1 )
499       ENDIF
500
501!
502!--    Change the droplet radius
[1359]503       IF ( ( new_r(n) - particles(n)%radius ) < 0.0_wp  .AND.                 &
504            new_r(n) < 0.0_wp )  THEN
505          WRITE( message_string, * ) '#1 k=',k,' j=',j,' i=',i,                &
[1849]506                       ' e_s=',e_s, ' e_a=',e_a,' t_int=',t_int,      &
[1359]507                       ' &delta_r=',delta_r,                                   &
[849]508                       ' particle_radius=',particles(n)%radius
509          CALL message( 'lpm_droplet_condensation', 'PA0144', 2, 2, -1, 6, 1 )
510       ENDIF
511
512!
513!--    Sum up the total volume of liquid water (needed below for
514!--    re-calculating the weighting factors)
[1359]515       ql_v(k,j,i) = ql_v(k,j,i) + particles(n)%weight_factor * new_r(n)**3
[849]516
[1359]517       particles(n)%radius = new_r(n)
[849]518
519!
520!--    Determine radius class of the particle needed for collision
[1359]521       IF ( ( hall_kernel  .OR.  wang_kernel )  .AND.  use_kernel_tables )     &
[849]522       THEN
[1359]523          particles(n)%class = ( LOG( new_r(n) ) - rclass_lbound ) /           &
524                               ( rclass_ubound - rclass_lbound ) *             &
[849]525                               radius_classes
526          particles(n)%class = MIN( particles(n)%class, radius_classes )
527          particles(n)%class = MAX( particles(n)%class, 1 )
528       ENDIF
529
530    ENDDO
531
532    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'stop' )
533
534
535 END SUBROUTINE lpm_droplet_condensation
Note: See TracBrowser for help on using the repository browser.