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

Last change on this file since 1856 was 1852, checked in by hoffmann, 8 years ago

last commit documented

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