source: palm/trunk/SOURCE/lpm_droplet_collision.f90 @ 3183

Last change on this file since 3183 was 3039, checked in by schwenkel, 6 years ago

bugfix for lcm with grid stretching

  • Property svn:keywords set to Id
File size: 15.4 KB
RevLine 
[1682]1!> @file lpm_droplet_collision.f90
[2000]2!------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]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!
[2718]17! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[849]20! Current revisions:
21! ------------------
[2375]22!
23!
[1321]24! Former revisions:
25! -----------------
26! $Id: lpm_droplet_collision.f90 3039 2018-05-24 13:13:11Z suehring $
[3039]27! bugfix for lcm with grid stretching
28!
29! 2718 2018-01-02 08:49:38Z maronga
[2716]30! Corrected "Former revisions" section
31!
32! 2696 2017-12-14 17:12:51Z kanani
33! Change in file header (GPL part)
34!
35! 2375 2017-08-29 14:10:28Z schwenkel
[2375]36! Changed ONLY-dependencies
37!
38! 2312 2017-07-14 20:26:51Z hoffmann
[2312]39! Consideration of aerosol mass during collision. Average impact algorithm has
40! been removed.
41!
42! 2274 2017-06-09 13:27:48Z Giersch
[2274]43! Changed error messages
[1321]44!
[2274]45! 2123 2017-01-18 12:34:59Z hoffmann
46!
[2123]47! 2122 2017-01-18 12:22:54Z hoffmann
48! Some reformatting of the code.
49!
[2001]50! 2000 2016-08-20 18:09:15Z knoop
51! Forced header and separation lines into 80 columns
52!
[1885]53! 1884 2016-04-21 11:11:40Z hoffmann
54! Conservation of mass should only be checked if collisions took place.
55!
[1861]56! 1860 2016-04-13 13:21:28Z hoffmann
[2312]57! Interpolation of dissipation rate adjusted to more reasonable values.
[1861]58!
[1823]59! 1822 2016-04-07 07:49:42Z hoffmann
60! Integration of a new collision algortithm based on Shima et al. (2009) and
61! Soelch and Kaercher (2010) called all_or_nothing. The previous implemented
62! collision algorithm is called average_impact. Moreover, both algorithms are
63! now positive definit due to their construction, i.e., no negative weighting
64! factors should occur.
65!
[1683]66! 1682 2015-10-07 23:56:08Z knoop
[2312]67! Code annotations made doxygen readable
68!
[1360]69! 1359 2014-04-11 17:15:14Z hoffmann
[2312]70! New particle structure integrated.
[1360]71! Kind definition added to all floating point numbers.
72!
[1323]73! 1322 2014-03-20 16:38:49Z raasch
74! REAL constants defined as wp_kind
75!
[1321]76! 1320 2014-03-20 08:40:49Z raasch
[1320]77! ONLY-attribute added to USE-statements,
78! kind-parameters added to all INTEGER and REAL declaration statements,
79! kinds are defined in new module kinds,
80! revision history before 2012 removed,
81! comment fields (!:) to be used for variable explanations added to
82! all variable declaration statements
[1072]83!
[1093]84! 1092 2013-02-02 11:24:22Z raasch
85! unused variables removed
86!
[1072]87! 1071 2012-11-29 16:54:55Z franke
[1071]88! Calculation of Hall and Wang kernel now uses collision-coalescence formulation
89! proposed by Wang instead of the continuous collection equation (for more
90! information about new method see PALM documentation)
91! Bugfix: message identifiers added
[849]92!
[1037]93! 1036 2012-10-22 13:43:42Z raasch
94! code put under GPL (PALM 3.9)
95!
[850]96! 849 2012-03-15 10:35:09Z raasch
97! initial revision (former part of advec_particles)
[849]98!
[850]99!
[849]100! Description:
101! ------------
[1682]102!> Calculates change in droplet radius by collision. Droplet collision is
103!> calculated for each grid box seperately. Collision is parameterized by
[1822]104!> using collision kernels. Two different kernels are available:
[1682]105!> Hall kernel: Kernel from Hall (1980, J. Atmos. Sci., 2486-2507), which
106!>              considers collision due to pure gravitational effects.
107!> Wang kernel: Beside gravitational effects (treated with the Hall-kernel) also
108!>              the effects of turbulence on the collision are considered using
109!>              parameterizations of Ayala et al. (2008, New J. Phys., 10,
110!>              075015) and Wang and Grabowski (2009, Atmos. Sci. Lett., 10,
111!>              1-8). This kernel includes three possible effects of turbulence:
112!>              the modification of the relative velocity between the droplets,
113!>              the effect of preferential concentration, and the enhancement of
[2312]114!>              collision efficiencies.
[849]115!------------------------------------------------------------------------------!
[1682]116 SUBROUTINE lpm_droplet_collision (i,j,k)
[849]117
[1320]118    USE arrays_3d,                                                             &
[3039]119        ONLY:  diss, dzw, ql_v, ql_vp
[849]120
[1320]121    USE cloud_parameters,                                                      &
[2375]122        ONLY:  rho_l, rho_s
[1320]123
124    USE constants,                                                             &
125        ONLY:  pi
126
127    USE control_parameters,                                                    &
[1822]128        ONLY:  dt_3d, message_string, dz
[1320]129
130    USE cpulog,                                                                &
131        ONLY:  cpu_log, log_point_s
132
133    USE grid_variables,                                                        &
[1822]134        ONLY:  dx, dy
[1320]135
136    USE kinds
137
138    USE lpm_collision_kernels_mod,                                             &
[1822]139        ONLY:  ckernel, recalculate_kernel
[1320]140
141    USE particle_attributes,                                                   &
[2312]142        ONLY:  curvature_solution_effects, dissipation_classes, hall_kernel,   &
143               iran_part, number_of_particles, particles, particle_type,       &
[2375]144               prt_count, use_kernel_tables, wang_kernel
[1320]145
[1822]146    USE random_function_mod,                                                   &
147        ONLY:  random_function
148
[1359]149    USE pegrid
150
[849]151    IMPLICIT NONE
152
[1682]153    INTEGER(iwp) ::  eclass   !<
154    INTEGER(iwp) ::  i        !<
155    INTEGER(iwp) ::  j        !<
156    INTEGER(iwp) ::  k        !<
157    INTEGER(iwp) ::  n        !<
[1822]158    INTEGER(iwp) ::  m        !<
[1682]159    INTEGER(iwp) ::  rclass_l !<
160    INTEGER(iwp) ::  rclass_s !<
[849]161
[1822]162    REAL(wp) ::  collection_probability  !< probability for collection
163    REAL(wp) ::  ddV                     !< inverse grid box volume
164    REAL(wp) ::  epsilon                 !< dissipation rate
165    REAL(wp) ::  factor_volume_to_mass   !< 4.0 / 3.0 * pi * rho_l
[2312]166    REAL(wp) ::  xm                      !< droplet mass of super-droplet m
167    REAL(wp) ::  xn                      !< droplet mass of super-droplet n
168    REAL(wp) ::  xsm                     !< aerosol mass of super-droplet m
169    REAL(wp) ::  xsn                     !< aerosol mass of super-droplet n
[1359]170
[2312]171    REAL(wp), DIMENSION(:), ALLOCATABLE ::  weight    !< weighting factor
172    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mass      !< total mass of super droplet
173    REAL(wp), DIMENSION(:), ALLOCATABLE ::  aero_mass !< total aerosol mass of super droplet
[1359]174
[849]175    CALL cpu_log( log_point_s(43), 'lpm_droplet_coll', 'start' )
176
[1822]177    number_of_particles   = prt_count(k,j,i)
[2312]178    factor_volume_to_mass = 4.0_wp / 3.0_wp * pi * rho_l
[3039]179    ddV                   = 1.0_wp / ( dx * dy * dzw(k) )
[849]180!
[1822]181!-- Collision requires at least one super droplet inside the box
182    IF ( number_of_particles > 0 )  THEN
[849]183
[1822]184       IF ( use_kernel_tables )  THEN
[849]185!
[1822]186!--       Fast method with pre-calculated collection kernels for
[1359]187!--       discrete radius- and dissipation-classes.
188          IF ( wang_kernel )  THEN
[1860]189             eclass = INT( diss(k,j,i) * 1.0E4_wp / 600.0_wp * &
[1359]190                           dissipation_classes ) + 1
191             epsilon = diss(k,j,i)
192          ELSE
193             epsilon = 0.0_wp
194          ENDIF
[2312]195
[1359]196          IF ( hall_kernel  .OR.  epsilon * 1.0E4_wp < 0.001_wp )  THEN
197             eclass = 0   ! Hall kernel is used
198          ELSE
199             eclass = MIN( dissipation_classes, eclass )
200          ENDIF
201
[2312]202       ELSE
[849]203!
[2312]204!--       Collection kernels are re-calculated for every new
[1822]205!--       grid box. First, allocate memory for kernel table.
206!--       Third dimension is 1, because table is re-calculated for
207!--       every new dissipation value.
208          ALLOCATE( ckernel(1:number_of_particles,1:number_of_particles,1:1) )
209!
210!--       Now calculate collection kernel for this box. Note that
211!--       the kernel is based on the previous time step
212          CALL recalculate_kernel( i, j, k )
[2312]213
214       ENDIF
[1822]215!
[2312]216!--    Temporary fields for total mass of super-droplet, aerosol mass, and
217!--    weighting factor are allocated.
218       ALLOCATE(mass(1:number_of_particles), weight(1:number_of_particles))
219       IF ( curvature_solution_effects )  ALLOCATE(aero_mass(1:number_of_particles))
[849]220
[2312]221       mass(1:number_of_particles)   = particles(1:number_of_particles)%weight_factor * &
222                                       particles(1:number_of_particles)%radius**3     * &
223                                       factor_volume_to_mass
[849]224
[2312]225       weight(1:number_of_particles) = particles(1:number_of_particles)%weight_factor
[849]226
[2312]227       IF ( curvature_solution_effects )  THEN
228          aero_mass(1:number_of_particles) = particles(1:number_of_particles)%weight_factor * &
229                                             particles(1:number_of_particles)%aux1**3       * &
230                                             4.0 / 3.0 * pi * rho_s
231       ENDIF
232!
233!--    Calculate collision/coalescence
234       DO  n = 1, number_of_particles
[849]235
[2312]236          DO  m = n, number_of_particles
237!
238!--          For collisions, the weighting factor of at least one super-droplet
239!--          needs to be larger or equal to one.
240             IF ( MIN( weight(n), weight(m) ) .LT. 1.0 )  CYCLE
241!
242!--          Get mass of individual droplets (aerosols)
243             xn = mass(n) / weight(n)
244             xm = mass(m) / weight(m)
245             IF ( curvature_solution_effects )  THEN
246                xsn = aero_mass(n) / weight(n)
247                xsm = aero_mass(m) / weight(m)
248             ENDIF
249!
250!--          Probability that the necessary collisions take place
251             IF ( use_kernel_tables )  THEN
252                rclass_l = particles(n)%class
253                rclass_s = particles(m)%class
[849]254
[2312]255                collection_probability  = MAX( weight(n), weight(m) ) *     &
256                                          ckernel(rclass_l,rclass_s,eclass) * ddV * dt_3d
257             ELSE
258                collection_probability  = MAX( weight(n), weight(m) ) *     &
259                                          ckernel(n,m,1) * ddV * dt_3d
260             ENDIF
[849]261!
[2312]262!--          Calculate the number of collections and consider multiple collections.
263!--          (Accordingly, p_crit will be 0.0, 1.0, 2.0, ...)
264             IF ( collection_probability - FLOOR(collection_probability)    &
265                  .GT. random_function( iran_part ) )  THEN
266                collection_probability = FLOOR(collection_probability) + 1.0_wp
267             ELSE
268                collection_probability = FLOOR(collection_probability)
269             ENDIF
[1822]270
[2312]271             IF ( collection_probability .GT. 0.0 )  THEN
[849]272!
[2312]273!--             Super-droplet n collects droplets of super-droplet m
274                IF ( weight(n) .LT. weight(m) )  THEN
[849]275
[2312]276                   mass(n)   = mass(n)   + weight(n) * xm * collection_probability
277                   weight(m) = weight(m) - weight(n)      * collection_probability
278                   mass(m)   = mass(m)   - weight(n) * xm * collection_probability
279                   IF ( curvature_solution_effects )  THEN
280                      aero_mass(n) = aero_mass(n) + weight(n) * xsm * collection_probability
281                      aero_mass(m) = aero_mass(m) - weight(n) * xsm * collection_probability
[1822]282                   ENDIF
[849]283
[2312]284                ELSEIF ( weight(m) .LT. weight(n) )  THEN
[849]285
[2312]286                   mass(m)   = mass(m)   + weight(m) * xn * collection_probability
287                   weight(n) = weight(n) - weight(m)      * collection_probability
288                   mass(n)   = mass(n)   - weight(m) * xn * collection_probability
289                   IF ( curvature_solution_effects )  THEN
290                      aero_mass(m) = aero_mass(m) + weight(m) * xsn * collection_probability
291                      aero_mass(n) = aero_mass(n) - weight(m) * xsn * collection_probability
292                   ENDIF
[849]293
[2312]294                ELSE
[1822]295!
[2312]296!--                Collisions of particles of the same weighting factor.
297!--                Particle n collects 1/2 weight(n) droplets of particle m,
298!--                particle m collects 1/2 weight(m) droplets of particle n.
299!--                The total mass mass changes accordingly.
300!--                If n = m, the first half of the droplets coalesces with the
301!--                second half of the droplets; mass is unchanged because
302!--                xm = xn for n = m.
[1822]303!--
[2312]304!--                Note: For m = n this equation is an approximation only
305!--                valid for weight >> 1 (which is usually the case). The
306!--                approximation is weight(n)-1 = weight(n).
307                   mass(n)   = mass(n)   + 0.5_wp * weight(n) * ( xm - xn )
308                   mass(m)   = mass(m)   + 0.5_wp * weight(m) * ( xn - xm )
309                   IF ( curvature_solution_effects )  THEN
310                      aero_mass(n) = aero_mass(n) + 0.5_wp * weight(n) * ( xsm - xsn )
311                      aero_mass(m) = aero_mass(m) + 0.5_wp * weight(m) * ( xsn - xsm )
[1822]312                   ENDIF
[2312]313                   weight(n) = weight(n) - 0.5_wp * weight(m)
314                   weight(m) = weight(n)
[1071]315
[2312]316                ENDIF
[1071]317
[2312]318             ENDIF
[1071]319
[2312]320          ENDDO
[1071]321
[2312]322          ql_vp(k,j,i) = ql_vp(k,j,i) + mass(n) / factor_volume_to_mass
[1071]323
[2312]324       ENDDO
[1071]325
[2312]326       IF ( ANY(weight < 0.0_wp) )  THEN
327             WRITE( message_string, * ) 'negative weighting factor'
328             CALL message( 'lpm_droplet_collision', 'PA0028',      &
329                            2, 2, -1, 6, 1 )
330       ENDIF
[1822]331
[2312]332       particles(1:number_of_particles)%radius = ( mass(1:number_of_particles) /   &
333                                                   ( weight(1:number_of_particles) &
334                                                     * factor_volume_to_mass       &
335                                                   )                               &
336                                                 )**0.33333333333333_wp
[1822]337
[2312]338       IF ( curvature_solution_effects )  THEN
339          particles(1:number_of_particles)%aux1 = ( aero_mass(1:number_of_particles) / &
340                                                    ( weight(1:number_of_particles)    &
341                                                      * 4.0_wp / 3.0_wp * pi * rho_s   &
342                                                    )                                  &
343                                                  )**0.33333333333333_wp
344       ENDIF
[1822]345
[2312]346       particles(1:number_of_particles)%weight_factor = weight(1:number_of_particles)
[1822]347
[2312]348       DEALLOCATE( weight, mass )
349       IF ( curvature_solution_effects )  DEALLOCATE( aero_mass )
350       IF ( .NOT. use_kernel_tables )  DEALLOCATE( ckernel )
351
[849]352!
[1884]353!--    Check if LWC is conserved during collision process
354       IF ( ql_v(k,j,i) /= 0.0_wp )  THEN
355          IF ( ql_vp(k,j,i) / ql_v(k,j,i) >= 1.0001_wp  .OR.                      &
356               ql_vp(k,j,i) / ql_v(k,j,i) <= 0.9999_wp )  THEN
357             WRITE( message_string, * ) ' LWC is not conserved during',           &
358                                        ' collision! ',                           &
359                                        ' LWC after condensation: ', ql_v(k,j,i), &
360                                        ' LWC after collision: ', ql_vp(k,j,i)
361             CALL message( 'lpm_droplet_collision', 'PA0040', 2, 2, -1, 6, 1 )
362          ENDIF
[1359]363       ENDIF
[1884]364
[1359]365    ENDIF
[2312]366
[849]367    CALL cpu_log( log_point_s(43), 'lpm_droplet_coll', 'stop' )
368
[2312]369 END SUBROUTINE lpm_droplet_collision
Note: See TracBrowser for help on using the repository browser.