source: palm/trunk/SOURCE/lpm_collision_kernels.f90 @ 849

Last change on this file since 849 was 849, checked in by raasch, 12 years ago

Changed:


Original routine advec_particles split into several new subroutines and renamed
lpm.
init_particles renamed lpm_init
user_advec_particles renamed user_lpm_advec,
particle_boundary_conds renamed lpm_boundary_conds,
set_particle_attributes renamed lpm_set_attributes,
user_init_particles renamed user_lpm_init,
user_particle_attributes renamed user_lpm_set_attributes
(Makefile, lpm_droplet_collision, lpm_droplet_condensation, init_3d_model, modules, palm, read_var_list, time_integration, write_var_list, deleted: advec_particles, init_particles, particle_boundary_conds, set_particle_attributes, user_advec_particles, user_init_particles, user_particle_attributes, new: lpm, lpm_advec, lpm_boundary_conds, lpm_calc_liquid_water_content, lpm_data_output_particles, lpm_droplet_collision, lpm_drollet_condensation, lpm_exchange_horiz, lpm_extend_particle_array, lpm_extend_tails, lpm_extend_tail_array, lpm_init, lpm_init_sgs_tke, lpm_pack_arrays, lpm_read_restart_file, lpm_release_set, lpm_set_attributes, lpm_sort_arrays, lpm_write_exchange_statistics, lpm_write_restart_file, user_lpm_advec, user_lpm_init, user_lpm_set_attributes

  • Property svn:keywords set to Id
File size: 35.3 KB
RevLine 
[828]1 MODULE lpm_collision_kernels_mod
[790]2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
[849]6! routine collision_efficiency_rogers added (moved from former advec_particles
7! to here)
[829]8!
9! Former revisions:
10! -----------------
11! $Id: lpm_collision_kernels.f90 849 2012-03-15 10:35:09Z raasch $
12!
[836]13! 835 2012-02-22 11:21:19Z raasch $
14! Bugfix: array diss can be used only in case of Wang kernel
15!
[829]16! 828 2012-02-21 12:00:36Z raasch
[828]17! code has been completely reformatted, routine colker renamed
18! recalculate_kernel,
19! routine init_kernels added, radius is now communicated to the collision
20! routines by array radclass
[790]21!
[828]22! Bugfix: transformation factor for dissipation changed from 1E5 to 1E4
23!
[826]24! 825 2012-02-19 03:03:44Z raasch
25! routine renamed from wang_kernel to lpm_collision_kernels,
26! turbulence_effects on collision replaced by wang_kernel
27!
[800]28! 799 2011-12-21 17:48:03Z franke
29! speed optimizations and formatting
30! Bugfix: iq=1 is not allowed (routine effic)
31! Bugfix: replaced stop by ec=0.0 in case of very small ec (routine effic)
32!
[791]33! 790 2011-11-29 03:11:20Z raasch
34! initial revision
[790]35!
36! Description:
37! ------------
[828]38! This module calculates collision efficiencies either due to pure gravitational
39! effects (Hall kernel, see Hall, 1980: J. Atmos. Sci., 2486-2507) or
40! including the effects of (SGS) turbulence (Wang kernel, see Wang and
41! Grabowski, 2009: Atmos. Sci. Lett., 10, 1-8). The original code has been
42! provided by L.-P. Wang but is substantially reformatted and speed optimized
43! here.
44!
45! ATTENTION:
46! Physical quantities (like g, densities, etc.) used in this module still
47! have to be adjusted to those values used in the main PALM code.
48! Also, quantities in CGS-units should be converted to SI-units eventually.
[790]49!------------------------------------------------------------------------------!
50
51    USE arrays_3d
52    USE cloud_parameters
53    USE constants
54    USE particle_attributes
[828]55    USE pegrid
[790]56
[828]57
[790]58    IMPLICIT NONE
59
60    PRIVATE
61
[849]62    PUBLIC  ckernel, collision_efficiency_rogers, init_kernels, &
63            rclass_lbound, rclass_ubound, recalculate_kernel           
[790]64
[828]65    REAL ::  epsilon, eps2, rclass_lbound, rclass_ubound, urms, urms2
[790]66
[828]67    REAL, DIMENSION(:),   ALLOCATABLE ::  epsclass, radclass, winf
68    REAL, DIMENSION(:,:), ALLOCATABLE ::  ec, ecf, gck, hkernel, hwratio
69    REAL, DIMENSION(:,:,:), ALLOCATABLE ::  ckernel
[792]70
[828]71    SAVE
[792]72
[790]73!
74!-- Public interfaces
[849]75    INTERFACE collision_efficiency_rogers
76       MODULE PROCEDURE collision_efficiency_rogers
77    END INTERFACE collision_efficiency_rogers
78
[828]79    INTERFACE init_kernels
80       MODULE PROCEDURE init_kernels
81    END INTERFACE init_kernels
[790]82
[828]83    INTERFACE recalculate_kernel
84       MODULE PROCEDURE recalculate_kernel
85    END INTERFACE recalculate_kernel
[790]86
87
[828]88    CONTAINS
[790]89
[792]90
[828]91    SUBROUTINE init_kernels
92!------------------------------------------------------------------------------!
93! Initialization of the collision efficiency matrix with fixed radius and
94! dissipation classes, calculated at simulation start only.
95!------------------------------------------------------------------------------!
[792]96
[828]97       IMPLICIT NONE
[792]98
[828]99       INTEGER ::  i, j, k
[790]100
[828]101
102!
103!--    Calculate collision efficiencies for fixed radius- and dissipation
104!--    classes
105       IF ( collision_kernel(6:9) == 'fast' )  THEN
106
107          ALLOCATE( ckernel(1:radius_classes,1:radius_classes,               &
108                    0:dissipation_classes), epsclass(1:dissipation_classes), &
109                    radclass(1:radius_classes) )
110
111!
112!--       Calculate the radius class bounds with logarithmic distances
113!--       in the interval [1.0E-6, 2.0E-4] m
114          rclass_lbound = LOG( 1.0E-6 )
115          rclass_ubound = LOG( 2.0E-4 )
116          radclass(1)   = 1.0E-6
117          DO  i = 2, radius_classes
118             radclass(i) = EXP( rclass_lbound +                                &
119                                ( rclass_ubound - rclass_lbound ) * ( i-1.0 ) /&
120                                ( radius_classes - 1.0 ) )
121!             IF ( myid == 0 )  THEN
122!                PRINT*, 'i=', i, ' r = ', radclass(i)*1.0E6
123!             ENDIF
124          ENDDO
125!
126!--       Collision routines expect radius to be in cm
127          radclass = radclass * 100.0
128
129!
130!--       Set the class bounds for dissipation in interval [0, 1000] cm**2/s**3
131          DO  i = 1, dissipation_classes
132             epsclass(i) = 1000.0 * REAL( i ) / dissipation_classes
133!             IF ( myid == 0 )  THEN
134!                PRINT*, 'i=', i, ' eps = ', epsclass(i)
135!             ENDIF
136          ENDDO
137!
138!--       Calculate collision efficiencies of the Wang/ayala kernel
139          ALLOCATE( ec(1:radius_classes,1:radius_classes),  &
140                    ecf(1:radius_classes,1:radius_classes), &
141                    gck(1:radius_classes,1:radius_classes), &
142                    winf(1:radius_classes) )
143
144          DO  k = 1, dissipation_classes
145
146             epsilon = epsclass(k)
147             urms    = 202.0 * ( epsilon / 400.0 )**( 1.0 / 3.0 )
148
149             CALL turbsd
150             CALL turb_enhance_eff
151             CALL effic
152
153             DO  j = 1, radius_classes
154                DO  i = 1, radius_classes
155                   ckernel(i,j,k) = ec(i,j) * gck(i,j) * ecf(i,j)
156                ENDDO
157             ENDDO
158
159          ENDDO
160
161!
162!--       Calculate collision efficiencies of the Hall kernel
163          ALLOCATE( hkernel(1:radius_classes,1:radius_classes), &
164                    hwratio(1:radius_classes,1:radius_classes) )
165
166          CALL fallg
167          CALL effic
168
169          DO  j = 1, radius_classes
170             DO  i =  1, radius_classes
171                hkernel(i,j) = pi * ( radclass(j) + radclass(i) )**2 &
172                                  * ec(i,j) * ABS( winf(j) - winf(i) )
173                ckernel(i,j,0) = hkernel(i,j)  ! hall kernel stored on index 0
174              ENDDO
175          ENDDO
176
177!
178!--       Test output of efficiencies
179          IF ( j == -1 )  THEN
180
181             PRINT*, '*** Hall kernel'
182             WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i)*1.0E4, i = 1,radius_classes )
183             DO  j = 1, radius_classes
184                WRITE ( *,'(F4.0,1X,20(F4.2,1X))' ) radclass(j), ( hkernel(i,j), i = 1,radius_classes )
185             ENDDO
186
187             DO  k = 1, dissipation_classes
188                DO  i = 1, radius_classes
189                   DO  j = 1, radius_classes
190                      IF ( hkernel(i,j) == 0.0 )  THEN
191                         hwratio(i,j) = 9999999.9
192                      ELSE
193                         hwratio(i,j) = ckernel(i,j,k) / hkernel(i,j)
194                      ENDIF
195                   ENDDO
196                ENDDO
197
198                PRINT*, '*** epsilon = ', epsclass(k)
199                WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i)*1.0E4, i = 1,radius_classes )
200                DO  j = 1, radius_classes
201!                   WRITE ( *,'(F4.0,1X,20(F4.2,1X))' ) radclass(j)*1.0E4, ( ckernel(i,j,k), i = 1,radius_classes )
202                   WRITE ( *,'(F4.0,1X,20(F4.2,1X))' ) radclass(j)*1.0E4, ( hwratio(i,j), i = 1,radius_classes )
203                ENDDO
204             ENDDO
205
206          ENDIF
207
208          DEALLOCATE( ec, ecf, epsclass, gck, hkernel, winf )
209
210          ckernel = ckernel * 1.0E-6   ! kernel is needed in m**3/s
211
212       ELSEIF( collision_kernel == 'hall'  .OR.  collision_kernel == 'wang' ) &
213       THEN
214!
215!--       Initial settings for Hall- and Wang-Kernel
216!--       To be done: move here parts from turbsd, fallg, ecoll, etc.
217       ENDIF
218
219    END SUBROUTINE init_kernels
220
221
[790]222!------------------------------------------------------------------------------!
[828]223! Calculation of collision kernels during each timestep and for each grid box
[790]224!------------------------------------------------------------------------------!
[828]225    SUBROUTINE recalculate_kernel( i1, j1, k1 )
[790]226
227       USE arrays_3d
228       USE cloud_parameters
229       USE constants
[792]230       USE cpulog
[790]231       USE indices
[792]232       USE interfaces
[790]233       USE particle_attributes
234
235       IMPLICIT NONE
236
[828]237       INTEGER ::  i, i1, j, j1, k1, pend, pstart
[790]238
239
[828]240       pstart = prt_start_index(k1,j1,i1)
241       pend   = prt_start_index(k1,j1,i1) + prt_count(k1,j1,i1) - 1
242       radius_classes = prt_count(k1,j1,i1)
[792]243
[828]244       ALLOCATE( ec(1:radius_classes,1:radius_classes), &
245                 radclass(1:radius_classes), winf(1:radius_classes) )
[790]246
[828]247!
248!--    Store particle radii on the radclass array. Collision routines
249!--    expect radii to be in cm.
250       radclass(1:radius_classes) = particles(pstart:pend)%radius * 100.0
[790]251
[835]252       IF ( wang_kernel )  THEN
253          epsilon = diss(k1,j1,i1) * 1.0E4   ! dissipation rate in cm**2/s**-3
254       ELSE
255          epsilon = 0.0
256       ENDIF
[828]257       urms    = 202.0 * ( epsilon / 400.0 )**( 0.33333333333 )
[790]258
[828]259       IF ( wang_kernel  .AND.  epsilon > 0.001 )  THEN
260!
261!--       Call routines to calculate efficiencies for the Wang kernel
262          ALLOCATE( gck(1:radius_classes,1:radius_classes), &
263                    ecf(1:radius_classes,1:radius_classes) )
[790]264
[828]265          CALL turbsd
266          CALL turb_enhance_eff
267          CALL effic
[790]268
[828]269          DO  j = 1, radius_classes
270             DO  i =  1, radius_classes
271                ckernel(pstart+i-1,pstart+j-1,1) = ec(i,j) * gck(i,j) * ecf(i,j)
[790]272             ENDDO
[828]273          ENDDO
[790]274
[828]275          DEALLOCATE( gck, ecf )
[790]276
277       ELSE
[828]278!
279!--       Call routines to calculate efficiencies for the Hall kernel
[790]280          CALL fallg
281          CALL effic
282
[828]283          DO  j = 1, radius_classes
284             DO  i =  1, radius_classes
285                ckernel(pstart+i-1,pstart+j-1,1) = pi *                       &
286                                          ( radclass(j) + radclass(i) )**2    &
287                                          * ec(i,j) * ABS( winf(j) - winf(i) )
[790]288             ENDDO
289          ENDDO
290
291       ENDIF
292
[828]293       ckernel = ckernel * 1.0E-6   ! kernel is needed in m**3/s
[790]294
[828]295       DEALLOCATE( ec, radclass, winf )
[790]296
[828]297    END SUBROUTINE recalculate_kernel
[790]298
[828]299
[790]300!------------------------------------------------------------------------------!
[828]301! Calculation of gck
302! This is from Aayala 2008b, page 37ff.
303! Necessary input parameters: water density, radii of droplets, air density,
304! air viscosity, turbulent dissipation rate, taylor microscale reynolds number,
305! gravitational acceleration  --> to be replaced by PALM parameters
[790]306!------------------------------------------------------------------------------!
[792]307    SUBROUTINE turbsd
[799]308
[790]309       USE constants
310       USE cloud_parameters
311       USE particle_attributes
312       USE arrays_3d
313
314       IMPLICIT NONE
315
[828]316       INTEGER ::  i, j
[790]317
[828]318       LOGICAL, SAVE ::  first = .TRUE.
[790]319
[828]320       REAL ::  ao, ao_gr, bbb, be, b1, b2, ccc, c1, c1_gr, c2, d1, d2, eta, &
321                e1, e2, fao_gr, fr, grfin, lambda, lambda_re, lf, rc, rrp,   &
322                sst, tauk, tl, t2, tt, t1, vk, vrms1xy, vrms2xy, v1, v1v2xy, &
323                v1xysq, v2, v2xysq, wrfin, wrgrav2, wrtur2xy, xx, yy, z
[799]324
[828]325       REAL, SAVE ::  airdens, airvisc, anu, gravity, waterdens
[799]326
[828]327       REAL, DIMENSION(1:radius_classes) ::  st, tau
[790]328
[828]329
[799]330!
[828]331!--    Initial assignment of constants
[799]332       IF ( first )  THEN
[790]333
[799]334          first = .FALSE.
[828]335          airvisc   = 0.1818     ! dynamic viscosity in mg/cm*s
336          airdens   = 1.2250     ! air density in mg/cm**3
337          waterdens = 1000.0     ! water density in mg/cm**3
338          gravity   = 980.6650   ! in cm/s**2
[799]339          anu = airvisc/airdens  ! kinetic viscosity in cm**2/s
[790]340
[799]341       ENDIF   
[790]342
[828]343       lambda    = urms * SQRT( 15.0 * anu / epsilon )     ! in cm
344       lambda_re = urms**2 * SQRT( 15.0 / epsilon / anu )
345       tl        = urms**2 / epsilon                       ! in s
346       lf        = 0.5 * urms**3 / epsilon                 ! in cm
347       tauk      = SQRT( anu / epsilon )                   ! in s
348       eta       = ( anu**3 / epsilon )**0.25              ! in cm
349       vk        = eta / tauk                              ! in cm/s
[790]350
[828]351       ao = ( 11.0 + 7.0 * lambda_re ) / ( 205.0 + lambda_re )
352       tt = SQRT( 2.0 * lambda_re / ( SQRT( 15.0 ) * ao ) ) * tauk   ! in s
[799]353
[828]354       CALL fallg    ! gives winf in cm/s
[790]355
[828]356       DO  i = 1, radius_classes
357          tau(i) = winf(i) / gravity    ! in s
358          st(i)  = tau(i) / tauk
[790]359       ENDDO
360
[828]361!
362!--    Calculate wr (from Aayala 2008b, page 38f)
363       z   = tt / tl
364       be  = SQRT( 2.0 ) * lambda / lf
365       bbb = SQRT( 1.0 - 2.0 * be**2 )
366       d1  = ( 1.0 + bbb ) / ( 2.0 * bbb )
367       e1  = lf * ( 1.0 + bbb ) * 0.5   ! in cm
368       d2  = ( 1.0 - bbb ) * 0.5 / bbb
369       e2  = lf * ( 1.0 - bbb ) * 0.5   ! in cm
370       ccc = SQRT( 1.0 - 2.0 * z**2 )
371       b1  = ( 1.0 + ccc ) * 0.5 / ccc
372       c1  = tl * ( 1.0 + ccc ) * 0.5   ! in s
373       b2  = ( 1.0 - ccc ) * 0.5 / ccc
374       c2  = tl * ( 1.0 - ccc ) * 0.5   ! in s
[790]375
[828]376       DO  i = 1, radius_classes
[790]377
[828]378          v1 = winf(i)        ! in cm/s
379          t1 = tau(i)         ! in s
[790]380
[828]381          DO  j = 1, i
382             rrp = radclass(i) + radclass(j)               ! radius in cm
383             v2  = winf(j)                                 ! in cm/s
384             t2  = tau(j)                                  ! in s
[790]385
[828]386             v1xysq  = b1 * d1 * phi(c1,e1,v1,t1) - b1 * d2 * phi(c1,e2,v1,t1) &
387                     - b2 * d1 * phi(c2,e1,v1,t1) + b2 * d2 * phi(c2,e2,v1,t1)
388             v1xysq  = v1xysq * urms**2 / t1                ! in cm**2/s**2
389             vrms1xy = SQRT( v1xysq )                       ! in cm/s
[790]390
[828]391             v2xysq  = b1 * d1 * phi(c1,e1,v2,t2) - b1 * d2 * phi(c1,e2,v2,t2) &
392                     - b2 * d1 * phi(c2,e1,v2,t2) + b2 * d2 * phi(c2,e2,v2,t2)
393             v2xysq  = v2xysq * urms**2 / t2                ! in cm**2/s**2
394             vrms2xy = SQRT( v2xysq )                       ! in cm/s
[790]395
[828]396             IF ( winf(i) >= winf(j) )  THEN
[799]397                v1 = winf(i)
[790]398                t1 = tau(i)
[799]399                v2 = winf(j)
[790]400                t2 = tau(j)
401             ELSE
[799]402                v1 = winf(j)
[790]403                t1 = tau(j)
[799]404                v2 = winf(i)
[790]405                t2 = tau(i)
406             ENDIF
407
[828]408             v1v2xy   =  b1 * d1 * zhi(c1,e1,v1,t1,v2,t2) - &
409                         b1 * d2 * zhi(c1,e2,v1,t1,v2,t2) - &
410                         b2 * d1 * zhi(c2,e1,v1,t1,v2,t2) + &
411                         b2 * d2* zhi(c2,e2,v1,t1,v2,t2)
412             fr       = d1 * EXP( -rrp / e1 ) - d2 * EXP( -rrp / e2 )
413             v1v2xy   = v1v2xy * fr * urms**2 / tau(i) / tau(j)  ! in cm**2/s**2
414             wrtur2xy = vrms1xy**2 + vrms2xy**2 - 2.0 * v1v2xy   ! in cm**2/s**2
415             IF ( wrtur2xy < 0.0 )  wrtur2xy = 0.0
416             wrgrav2  = pi / 8.0 * ( winf(j) - winf(i) )**2
417             wrfin    = SQRT( ( 2.0 / pi ) * ( wrtur2xy + wrgrav2) )   ! in cm/s
[790]418
[828]419!
420!--          Calculate gr
421             IF ( st(j) > st(i) )  THEN
422                sst = st(j)
[790]423             ELSE
[828]424                sst = st(i)
[790]425             ENDIF
426
[828]427             xx = -0.1988 * sst**4 + 1.5275 * sst**3 - 4.2942 * sst**2 + &
428                   5.3406 * sst
429             IF ( xx < 0.0 )  xx = 0.0
430             yy = 0.1886 * EXP( 20.306 / lambda_re )
[790]431
[828]432             c1_gr  =  xx / ( gravity / vk * tauk )**yy
[790]433
[828]434             ao_gr  = ao + ( pi / 8.0) * ( gravity / vk * tauk )**2
435             fao_gr = 20.115 * SQRT( ao_gr / lambda_re )
436             rc     = SQRT( fao_gr * ABS( st(j) - st(i) ) ) * eta   ! in cm
[790]437
[828]438             grfin  = ( ( eta**2 + rc**2 ) / ( rrp**2 + rc**2) )**( c1_gr*0.5 )
439             IF ( grfin < 1.0 )  grfin = 1.0
[790]440
[828]441             gck(i,j) = 2.0 * pi * rrp**2 * wrfin * grfin           ! in cm**3/s
[790]442             gck(j,i) = gck(i,j)
443
444          ENDDO
445       ENDDO
446
[828]447    END SUBROUTINE turbsd
[790]448
[828]449
[790]450!------------------------------------------------------------------------------!
[828]451! phi as a function
[790]452!------------------------------------------------------------------------------!
[828]453    REAL FUNCTION phi( a, b, vsett, tau0 )
[790]454
455       IMPLICIT NONE
456
[828]457       REAL ::  a, aa1, b, tau0, vsett
[790]458
[828]459       aa1 = 1.0 / tau0 + 1.0 / a + vsett / b
460       phi = 1.0 / aa1  - 0.5 * vsett / b / aa1**2  ! in s
[790]461
[828]462    END FUNCTION phi
[792]463
[790]464
465!------------------------------------------------------------------------------!
[828]466! zeta as a function
[790]467!------------------------------------------------------------------------------!
[828]468    REAL FUNCTION zhi( a, b, vsett1, tau1, vsett2, tau2 )
[790]469
470       IMPLICIT NONE
471
[828]472       REAL ::  a, aa1, aa2, aa3, aa4, aa5, aa6, b, tau1, tau2, vsett1, vsett2
[790]473
[828]474       aa1 = vsett2 / b - 1.0 / tau2 - 1.0 / a
475       aa2 = vsett1 / b + 1.0 / tau1 + 1.0 / a
476       aa3 = ( vsett1 - vsett2 ) / b + 1.0 / tau1 + 1.0 / tau2
477       aa4 = ( vsett2 / b )**2 - ( 1.0 / tau2 + 1.0 / a )**2
478       aa5 = vsett2 / b + 1.0 / tau2 + 1.0 / a
479       aa6 = 1.0 / tau1 - 1.0 / a + ( 1.0 / tau2 + 1.0 / a) * vsett1 / vsett2
480       zhi = (1.0 / aa1 - 1.0 / aa2 ) * ( vsett1 - vsett2 ) * 0.5 / b / aa3**2 &
481           + (4.0 / aa4 - 1.0 / aa5**2 - 1.0 / aa1**2 ) * vsett2 * 0.5 / b /aa6&
482           + (2.0 * ( b / aa2 - b / aa1 ) - vsett1 / aa2**2 + vsett2 / aa1**2 )&
483           * 0.5 / b / aa3      ! in s**2
[799]484
[828]485    END FUNCTION zhi
[790]486
[828]487
[790]488!------------------------------------------------------------------------------!
[828]489! Calculation of terminal velocity winf
[790]490!------------------------------------------------------------------------------!
[828]491    SUBROUTINE fallg
[790]492
493       USE constants
494       USE cloud_parameters
495       USE particle_attributes
496       USE arrays_3d
497
[828]498       IMPLICIT NONE
[790]499
[828]500       INTEGER ::  i, j
[790]501
[828]502       LOGICAL, SAVE ::  first = .TRUE.
[790]503
[828]504       REAL, SAVE ::  cunh, eta, grav, phy, py, rhoa, rhow, sigma, stb, stok, &
505                      t0, xlamb
[790]506
[828]507       REAL ::  bond, x, xrey, y
[799]508
[828]509       REAL, DIMENSION(1:7), SAVE  ::  b
510       REAL, DIMENSION(1:6), SAVE  ::  c
[799]511
512!
[828]513!--    Initial assignment of constants
514       IF ( first )  THEN
[799]515
[828]516          first = .FALSE.
517          b = (/  -0.318657E1,  0.992696E0, -0.153193E-2, -0.987059E-3, &
518                 -0.578878E-3, 0.855176E-4, -0.327815E-5 /)
519          c = (/  -0.500015E1,  0.523778E1,  -0.204914E1,   0.475294E0, &
520                 -0.542819E-1, 0.238449E-2 /)
[790]521
[828]522          eta   = 1.818E-4          ! in poise = g/(cm s)
523          xlamb = 6.62E-6           ! in cm
524          rhow  = 1.0               ! in g/cm**3
525          rhoa  = 1.225E-3          ! in g/cm**3
526          grav  = 980.665           ! in cm/s**2
527          cunh  = 1.257 * xlamb     ! in cm
528          t0    = 273.15            ! in K
529          sigma = 76.1 - 0.155 * ( 293.15 - t0 )              ! in N/m = g/s**2
530          stok  = 2.0  * grav * ( rhow - rhoa ) / ( 9.0 * eta ) ! in 1/(cm s)
531          stb   = 32.0 * rhoa * ( rhow - rhoa) * grav / (3.0 * eta * eta)
532                                                                ! in 1/cm**3
533          phy   = sigma**3 * rhoa**2 / ( eta**4 * grav * ( rhow - rhoa ) )
534          py    = phy**( 1.0 / 6.0 )
[790]535
[828]536       ENDIF
[790]537
[828]538       DO  j = 1, radius_classes
[790]539
[828]540          IF ( radclass(j) <= 1.0E-3 ) THEN
[799]541
[828]542             winf(j) = stok * ( radclass(j)**2 + cunh * radclass(j) ) ! in cm/s
[790]543
[828]544          ELSEIF ( radclass(j) > 1.0E-3  .AND.  radclass(j) <= 5.35E-2 )  THEN
[790]545
[828]546             x = LOG( stb * radclass(j)**3 )
547             y = 0.0
[790]548
[828]549             DO  i = 1, 7
550                y = y + b(i) * x**(i-1)
551             ENDDO
[790]552
[828]553             xrey = ( 1.0 + cunh / radclass(j) ) * EXP( y )
554             winf(j) = xrey * eta / ( 2.0 * rhoa * radclass(j) )      ! in cm/s
[790]555
[828]556          ELSEIF ( radclass(j) > 5.35E-2 )  THEN
[790]557
[828]558             IF ( radclass(j) > 0.35 )  THEN
559                bond = grav * ( rhow - rhoa ) * 0.35**2 / sigma
560             ELSE
561                bond = grav * ( rhow - rhoa ) * radclass(j)**2 / sigma
562             ENDIF
[790]563
[828]564             x = LOG( 16.0 * bond * py / 3.0 )
565             y = 0.0
[790]566
[828]567             DO  i = 1, 6
568                y = y + c(i) * x**(i-1)
569             ENDDO
[790]570
[828]571             xrey = py * EXP( y )
[790]572
[828]573             IF ( radclass(j) > 0.35 )  THEN
574                winf(j) = xrey * eta / ( 2.0 * rhoa * 0.35 )           ! in cm/s
575             ELSE
576                winf(j) = xrey * eta / ( 2.0 * rhoa * radclass(j) )    ! in cm/s
577             ENDIF
[790]578
[828]579          ENDIF
[790]580
[828]581       ENDDO
[790]582
[828]583    END SUBROUTINE fallg
[790]584
[828]585
[790]586!------------------------------------------------------------------------------!
[828]587! Calculation of collision efficencies for the Hall kernel
[790]588!------------------------------------------------------------------------------!
[828]589    SUBROUTINE effic
[790]590
[828]591       USE arrays_3d
592       USE cloud_parameters
593       USE constants
594       USE particle_attributes
[790]595
[828]596       IMPLICIT NONE
[790]597
[828]598       INTEGER ::  i, iq, ir, j, k, kk
[790]599
[828]600       INTEGER, DIMENSION(:), ALLOCATABLE ::  ira
[790]601
[828]602       LOGICAL, SAVE ::  first = .TRUE.
[790]603
[828]604       REAL ::  ek, particle_radius, pp, qq, rq
[790]605
[828]606       REAL, DIMENSION(1:21), SAVE ::  rat
607       REAL, DIMENSION(1:15), SAVE ::  r0
608       REAL, DIMENSION(1:15,1:21), SAVE ::  ecoll
[790]609
[792]610!
[828]611!--    Initial assignment of constants
612       IF ( first )  THEN
[790]613
[792]614         first = .FALSE.
[828]615         r0  = (/ 6.0, 8.0, 10.0, 15.0, 20.0, 25.0, 30.0, 40.0, 50.0, 60., &
616                  70.0, 100.0, 150.0, 200.0, 300.0 /)
617         rat = (/ 0.00, 0.05, 0.10, 0.15, 0.20, 0.25, 0.30, 0.35, 0.40, 0.45, &
618                  0.50, 0.55, 0.60, 0.65, 0.70, 0.75, 0.80, 0.85, 0.90, 0.95, &
619                  1.00 /)
[790]620
[828]621         ecoll(:,1) = (/0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, &
622                        0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001/)
623         ecoll(:,2) = (/0.003, 0.003, 0.003, 0.004, 0.005, 0.005, 0.005, &
624                        0.010, 0.100, 0.050, 0.200, 0.500, 0.770, 0.870, 0.970/)
625         ecoll(:,3) = (/0.007, 0.007, 0.007, 0.008, 0.009, 0.010, 0.010, &
626                        0.070, 0.400, 0.430, 0.580, 0.790, 0.930, 0.960, 1.000/)
627         ecoll(:,4) = (/0.009, 0.009, 0.009, 0.012, 0.015, 0.010, 0.020, &
628                        0.280, 0.600, 0.640, 0.750, 0.910, 0.970, 0.980, 1.000/)
629         ecoll(:,5) = (/0.014, 0.014, 0.014, 0.015, 0.016, 0.030, 0.060, &
630                        0.500, 0.700, 0.770, 0.840, 0.950, 0.970, 1.000, 1.000/)
631         ecoll(:,6) = (/0.017, 0.017, 0.017, 0.020, 0.022, 0.060, 0.100, &
632                        0.620, 0.780, 0.840, 0.880, 0.950, 1.000, 1.000, 1.000/)
633         ecoll(:,7) = (/0.030, 0.030, 0.024, 0.022, 0.032, 0.062, 0.200, &
634                        0.680, 0.830, 0.870, 0.900, 0.950, 1.000, 1.000, 1.000/)
635         ecoll(:,8) = (/0.025, 0.025, 0.025, 0.036, 0.043, 0.130, 0.270, &
636                        0.740, 0.860, 0.890, 0.920, 1.000, 1.000, 1.000, 1.000/)
637         ecoll(:,9) = (/0.027, 0.027, 0.027, 0.040, 0.052, 0.200, 0.400, &
638                        0.780, 0.880, 0.900, 0.940, 1.000, 1.000, 1.000, 1.000/)
639         ecoll(:,10)= (/0.030, 0.030, 0.030, 0.047, 0.064, 0.250, 0.500, &
640                        0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/)
641         ecoll(:,11)= (/0.040, 0.040, 0.033, 0.037, 0.068, 0.240, 0.550, &
642                        0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/)
643         ecoll(:,12)= (/0.035, 0.035, 0.035, 0.055, 0.079, 0.290, 0.580, &
644                        0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/)
645         ecoll(:,13)= (/0.037, 0.037, 0.037, 0.062, 0.082, 0.290, 0.590, &
646                        0.780, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/)
647         ecoll(:,14)= (/0.037, 0.037, 0.037, 0.060, 0.080, 0.290, 0.580, &
648                        0.770, 0.890, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/)
649         ecoll(:,15)= (/0.037, 0.037, 0.037, 0.041, 0.075, 0.250, 0.540, &
650                        0.760, 0.880, 0.920, 0.950, 1.000, 1.000, 1.000, 1.000/)
651         ecoll(:,16)= (/0.037, 0.037, 0.037, 0.052, 0.067, 0.250, 0.510, &
652                        0.770, 0.880, 0.930, 0.970, 1.000, 1.000, 1.000, 1.000/)
653         ecoll(:,17)= (/0.037, 0.037, 0.037, 0.047, 0.057, 0.250, 0.490, &
654                        0.770, 0.890, 0.950, 1.000, 1.000, 1.000, 1.000, 1.000/)
655         ecoll(:,18)= (/0.036, 0.036, 0.036, 0.042, 0.048, 0.230, 0.470, &
656                        0.780, 0.920, 1.000, 1.020, 1.020, 1.020, 1.020, 1.020/)
657         ecoll(:,19)= (/0.040, 0.040, 0.035, 0.033, 0.040, 0.112, 0.450, &
658                        0.790, 1.010, 1.030, 1.040, 1.040, 1.040, 1.040, 1.040/)
659         ecoll(:,20)= (/0.033, 0.033, 0.033, 0.033, 0.033, 0.119, 0.470, &
660                        0.950, 1.300, 1.700, 2.300, 2.300, 2.300, 2.300, 2.300/)
661         ecoll(:,21)= (/0.027, 0.027, 0.027, 0.027, 0.027, 0.125, 0.520, &
662                        1.400, 2.300, 3.000, 4.000, 4.000, 4.000, 4.000, 4.000/)
663       ENDIF
[790]664
[792]665!
[828]666!--    Calculate the radius class index of particles with respect to array r
667       ALLOCATE( ira(1:radius_classes) )
668       DO  j = 1, radius_classes
669          particle_radius = radclass(j) * 1.0E4   ! in microm
670          DO  k = 1, 15
671             IF ( particle_radius < r0(k) )  THEN
672                ira(j) = k
673                EXIT
674             ENDIF
675          ENDDO
676          IF ( particle_radius >= r0(15) )  ira(j) = 16
677       ENDDO
[790]678
[792]679!
[828]680!--    Two-dimensional linear interpolation of the collision efficiency.
681!--    Radius has to be in µm
682       DO  j = 1, radius_classes
683          DO  i = 1, j
[792]684
[828]685             ir = ira(j)
686             rq = radclass(i) / radclass(j)
687             iq = INT( rq * 20 ) + 1
688             iq = MAX( iq , 2)
[792]689
[828]690             IF ( ir < 16 )  THEN
691                IF ( ir >= 2 )  THEN
692                   pp = ( ( radclass(j) * 1.0E04 ) - r0(ir-1) ) / &
693                        ( r0(ir) - r0(ir-1) )
694                   qq = ( rq- rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
695                   ec(j,i) = ( 1.0-pp ) * ( 1.0-qq ) * ecoll(ir-1,iq-1)  &
696                             + pp * ( 1.0-qq ) * ecoll(ir,iq-1)          &
697                             + qq * ( 1.0-pp ) * ecoll(ir-1,iq)          &
698                             + pp * qq * ecoll(ir,iq)
699                ELSE
700                   qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
701                   ec(j,i) = (1.0-qq) * ecoll(1,iq-1) + qq * ecoll(1,iq)
702                ENDIF
703             ELSE
704                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
705                ek = ( 1.0 - qq ) * ecoll(15,iq-1) + qq * ecoll(15,iq)
706                ec(j,i) = MIN( ek, 1.0 )
707             ENDIF
[792]708
[828]709             ec(i,j) = ec(j,i)
710             IF ( ec(i,j) < 1.0E-20 )  ec(i,j) = 0.0
[792]711
[828]712          ENDDO
713       ENDDO
[792]714
[828]715       DEALLOCATE( ira )
[792]716
[828]717    END SUBROUTINE effic
[792]718
719
[790]720!------------------------------------------------------------------------------!
[828]721! Calculation of enhancement factor for collision efficencies due to turbulence
[790]722!------------------------------------------------------------------------------!
[828]723    SUBROUTINE turb_enhance_eff
[790]724
725       USE constants
726       USE cloud_parameters
727       USE particle_attributes
728       USE arrays_3d
729
[828]730       IMPLICIT NONE
[790]731
[828]732       INTEGER :: i, ik, iq, ir, j, k, kk
[790]733
[828]734       INTEGER, DIMENSION(:), ALLOCATABLE ::  ira
[790]735
[828]736       REAL ::  particle_radius, pp, qq, rq, x1, x2, x3, y1, y2, y3
[790]737
[828]738       LOGICAL, SAVE ::  first = .TRUE.
[799]739
[828]740       REAL, DIMENSION(1:11), SAVE ::  rat
741       REAL, DIMENSION(1:7), SAVE  ::  r0
742       REAL, DIMENSION(1:7,1:11), SAVE ::  ecoll_100, ecoll_400
[799]743
744!
[828]745!--    Initial assignment of constants
746       IF ( first )  THEN
[799]747
[828]748          first = .FALSE.
[799]749
[828]750          r0  = (/ 10.0, 20.0, 30.0, 40.0, 50.0, 60.0, 100.0 /)
751          rat = (/ 0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0 /)
752!
753!--       In 100 cm^2/s^3
754          ecoll_100(:,1) = (/1.74,  1.74,  1.773, 1.49,  1.207,  1.207,  1.0 /)
755          ecoll_100(:,2) = (/1.46,  1.46,  1.421, 1.245, 1.069,  1.069,  1.0 /)
756          ecoll_100(:,3) = (/1.32,  1.32,  1.245, 1.123, 1.000,  1.000,  1.0 /)
757          ecoll_100(:,4) = (/1.250, 1.250, 1.148, 1.087, 1.025,  1.025,  1.0 /)
758          ecoll_100(:,5) = (/1.186, 1.186, 1.066, 1.060, 1.056,  1.056,  1.0 /)
759          ecoll_100(:,6) = (/1.045, 1.045, 1.000, 1.014, 1.028,  1.028,  1.0 /)
760          ecoll_100(:,7) = (/1.070, 1.070, 1.030, 1.038, 1.046,  1.046,  1.0 /)
761          ecoll_100(:,8) = (/1.000, 1.000, 1.054, 1.042, 1.029,  1.029,  1.0 /)
762          ecoll_100(:,9) = (/1.223, 1.223, 1.117, 1.069, 1.021,  1.021,  1.0 /)
763          ecoll_100(:,10)= (/1.570, 1.570, 1.244, 1.166, 1.088,  1.088,  1.0 /)
764          ecoll_100(:,11)= (/20.3,  20.3,  14.6 , 8.61,  2.60,   2.60 ,  1.0 /)
765!
766!--       In 400 cm^2/s^3
767          ecoll_400(:,1) = (/4.976, 4.976,  3.593, 2.519, 1.445,  1.445,  1.0 /)
768          ecoll_400(:,2) = (/2.984, 2.984,  2.181, 1.691, 1.201,  1.201,  1.0 /)
769          ecoll_400(:,3) = (/1.988, 1.988,  1.475, 1.313, 1.150,  1.150,  1.0 /)
770          ecoll_400(:,4) = (/1.490, 1.490,  1.187, 1.156, 1.126,  1.126,  1.0 /)
771          ecoll_400(:,5) = (/1.249, 1.249,  1.088, 1.090, 1.092,  1.092,  1.0 /)
772          ecoll_400(:,6) = (/1.139, 1.139,  1.130, 1.091, 1.051,  1.051,  1.0 /)
773          ecoll_400(:,7) = (/1.220, 1.220,  1.190, 1.138, 1.086,  1.086,  1.0 /)
774          ecoll_400(:,8) = (/1.325, 1.325,  1.267, 1.165, 1.063,  1.063,  1.0 /)
775          ecoll_400(:,9) = (/1.716, 1.716,  1.345, 1.223, 1.100,  1.100,  1.0 /)
776          ecoll_400(:,10)= (/3.788, 3.788,  1.501, 1.311, 1.120,  1.120,  1.0 /)
777          ecoll_400(:,11)= (/36.52, 36.52,  19.16, 22.80,  26.0,   26.0,  1.0 /)
[799]778
[828]779       ENDIF
[790]780
[828]781!
782!--    Calculate the radius class index of particles with respect to array r0
783       ALLOCATE( ira(1:radius_classes) )
[790]784
[828]785       DO  j = 1, radius_classes
786          particle_radius = radclass(j) * 1.0E4    ! in microm
787          DO  k = 1, 7
788             IF ( particle_radius < r0(k) )  THEN
789                ira(j) = k
790                EXIT
791             ENDIF
792          ENDDO
793          IF ( particle_radius >= r0(7) )  ira(j) = 8
794       ENDDO
[799]795
796!
[828]797!--    Two-dimensional linear interpolation of the collision efficiencies
798       DO  j =  1, radius_classes
799          DO  i = 1, j
[799]800
[828]801             ir = ira(j)
802             rq = radclass(i) / radclass(j)
[799]803
[828]804             DO  kk = 2, 11
805                IF ( rq <= rat(kk) )  THEN
806                   iq = kk
807                   EXIT
808                ENDIF
809             ENDDO
[790]810
[828]811             y1 = 1.0      ! 0  cm2/s3
812!
813!--          100 cm2/s3, 400 cm2/s3
814             IF ( ir < 8 )  THEN
815                IF ( ir >= 2 )  THEN
816                   pp = ( radclass(j)*1.0E4 - r0(ir-1) ) / ( r0(ir) - r0(ir-1) )
817                   qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
818                   y2 = ( 1.0-pp ) * ( 1.0-qq ) * ecoll_100(ir-1,iq-1) +  &
819                                pp * ( 1.0-qq ) * ecoll_100(ir,iq-1)   +  &
820                                qq * ( 10.-pp ) * ecoll_100(ir-1,iq)   +  &
821                                pp * qq         * ecoll_100(ir,iq)
822                   y3 = ( 1.0-pp ) * ( 1.0-qq ) * ecoll_400(ir-1,iq-1) +  &
823                                pp * ( 1.0-qq ) * ecoll_400(ir,iq-1)   +  &
824                                qq * ( 1.0-pp ) * ecoll_400(ir-1,iq)   +  &
825                                pp * qq         * ecoll_400(ir,iq)
826                ELSE
827                   qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
828                   y2 = ( 1.0-qq ) * ecoll_100(1,iq-1) + qq * ecoll_100(1,iq)
829                   y3 = ( 1.0-qq ) * ecoll_400(1,iq-1) + qq * ecoll_400(1,iq)
830                ENDIF
831             ELSE
832                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
833                y2 = ( 1.0-qq ) * ecoll_100(7,iq-1) + qq * ecoll_100(7,iq)
834                y3 = ( 1.0-qq ) * ecoll_400(7,iq-1) + qq * ecoll_400(7,iq)
835             ENDIF
836!
837!--          Linear interpolation of dissipation rate in cm2/s3
838             IF ( epsilon <= 100.0 )  THEN
839                ecf(j,i) = ( epsilon - 100.0 ) / (   0.0 - 100.0 ) * y1 &
840                         + ( epsilon -   0.0 ) / ( 100.0 -   0.0 ) * y2
841             ELSEIF ( epsilon <= 600.0 )  THEN
842                ecf(j,i) = ( epsilon - 400.0 ) / ( 100.0 - 400.0 ) * y2 &
843                         + ( epsilon - 100.0 ) / ( 400.0 - 100.0 ) * y3
844             ELSE
845                ecf(j,i) = (   600.0 - 400.0 ) / ( 100.0 - 400.0 ) * y2 &
846                         + (   600.0 - 100.0 ) / ( 400.0 - 100.0 ) * y3
847             ENDIF
[790]848
[828]849             IF ( ecf(j,i) < 1.0 )  ecf(j,i) = 1.0
[790]850
[828]851             ecf(i,j) = ecf(j,i)
[790]852
[828]853          ENDDO
854       ENDDO
[790]855
[828]856    END SUBROUTINE turb_enhance_eff
[790]857
[849]858
859
860    SUBROUTINE collision_efficiency_rogers( mean_r, r, e)
861!------------------------------------------------------------------------------!
862! Collision efficiencies from table 8.2 in Rogers and Yau (1989, 3rd edition).
863! Values are calculated from table by bilinear interpolation.
864!------------------------------------------------------------------------------!
865
866       IMPLICIT NONE
867
868       INTEGER       ::  i, j, k
869
870       LOGICAL, SAVE ::  first = .TRUE.
871
872       REAL          ::  aa, bb, cc, dd, dx, dy, e, gg, mean_r, mean_rm, r, &
873                         rm, x, y
874
875       REAL, DIMENSION(1:9), SAVE      ::  collected_r = 0.0
876       REAL, DIMENSION(1:19), SAVE     ::  collector_r = 0.0
877       REAL, DIMENSION(1:9,1:19), SAVE ::  ef = 0.0
878
879       mean_rm = mean_r * 1.0E06
880       rm      = r      * 1.0E06
881
882       IF ( first )  THEN
883
884          collected_r = (/ 2.0, 3.0, 4.0, 6.0, 8.0, 10.0, 15.0, 20.0, 25.0 /)
885          collector_r = (/ 10.0, 20.0, 30.0, 40.0, 50.0, 60.0, 80.0, 100.0,  &
886                           150.0, 200.0, 300.0, 400.0, 500.0, 600.0, 1000.0, &
887                           1400.0, 1800.0, 2400.0, 3000.0 /)
888
889          ef(:,1) = (/0.017, 0.027, 0.037, 0.052, 0.052, 0.052, 0.052, 0.0, &
890                      0.0 /)
891          ef(:,2) = (/0.001, 0.016, 0.027, 0.060, 0.12, 0.17, 0.17, 0.17, 0.0 /)
892          ef(:,3) = (/0.001, 0.001, 0.02,  0.13,  0.28, 0.37, 0.54, 0.55, 0.47/)
893          ef(:,4) = (/0.001, 0.001, 0.02,  0.23,  0.4,  0.55, 0.7,  0.75, 0.75/)
894          ef(:,5) = (/0.01,  0.01,  0.03,  0.3,   0.4,  0.58, 0.73, 0.75, 0.79/)
895          ef(:,6) = (/0.01,  0.01,  0.13,  0.38,  0.57, 0.68, 0.80, 0.86, 0.91/)
896          ef(:,7) = (/0.01,  0.085, 0.23,  0.52,  0.68, 0.76, 0.86, 0.92, 0.95/)
897          ef(:,8) = (/0.01,  0.14,  0.32,  0.60,  0.73, 0.81, 0.90, 0.94, 0.96/)
898          ef(:,9) = (/0.025, 0.25,  0.43,  0.66,  0.78, 0.83, 0.92, 0.95, 0.96/)
899          ef(:,10)= (/0.039, 0.3,   0.46,  0.69,  0.81, 0.87, 0.93, 0.95, 0.96/)
900          ef(:,11)= (/0.095, 0.33,  0.51,  0.72,  0.82, 0.87, 0.93, 0.96, 0.97/)
901          ef(:,12)= (/0.098, 0.36,  0.51,  0.73,  0.83, 0.88, 0.93, 0.96, 0.97/)
902          ef(:,13)= (/0.1,   0.36,  0.52,  0.74,  0.83, 0.88, 0.93, 0.96, 0.97/)
903          ef(:,14)= (/0.17,  0.4,   0.54,  0.72,  0.83, 0.88, 0.94, 0.98, 1.0 /)
904          ef(:,15)= (/0.15,  0.37,  0.52,  0.74,  0.82, 0.88, 0.94, 0.98, 1.0 /)
905          ef(:,16)= (/0.11,  0.34,  0.49,  0.71,  0.83, 0.88, 0.94, 0.95, 1.0 /)
906          ef(:,17)= (/0.08,  0.29,  0.45,  0.68,  0.8,  0.86, 0.96, 0.94, 1.0 /)
907          ef(:,18)= (/0.04,  0.22,  0.39,  0.62,  0.75, 0.83, 0.92, 0.96, 1.0 /)
908          ef(:,19)= (/0.02,  0.16,  0.33,  0.55,  0.71, 0.81, 0.90, 0.94, 1.0 /)
909
910       ENDIF
911
912       DO  k = 1, 8
913          IF ( collected_r(k) <= mean_rm )  i = k
914       ENDDO
915
916       DO  k = 1, 18
917          IF ( collector_r(k) <= rm )  j = k
918       ENDDO
919
920       IF ( rm < 10.0 )  THEN
921          e = 0.0
922       ELSEIF ( mean_rm < 2.0 )  THEN
923          e = 0.001
924       ELSEIF ( mean_rm >= 25.0 )  THEN
925          IF( j <= 2 )  e = 0.0
926          IF( j == 3 )  e = 0.47
927          IF( j == 4 )  e = 0.8
928          IF( j == 5 )  e = 0.9
929          IF( j >=6  )  e = 1.0
930       ELSEIF ( rm >= 3000.0 )  THEN
931          IF( i == 1 )  e = 0.02
932          IF( i == 2 )  e = 0.16
933          IF( i == 3 )  e = 0.33
934          IF( i == 4 )  e = 0.55
935          IF( i == 5 )  e = 0.71
936          IF( i == 6 )  e = 0.81
937          IF( i == 7 )  e = 0.90
938          IF( i >= 8 )  e = 0.94
939       ELSE
940          x  = mean_rm - collected_r(i)
941          y  = rm - collector_r(j)
942          dx = collected_r(i+1) - collected_r(i)
943          dy = collector_r(j+1) - collector_r(j)
944          aa = x**2 + y**2
945          bb = ( dx - x )**2 + y**2
946          cc = x**2 + ( dy - y )**2
947          dd = ( dx - x )**2 + ( dy - y )**2
948          gg = aa + bb + cc + dd
949
950          e = ( (gg-aa)*ef(i,j) + (gg-bb)*ef(i+1,j) + (gg-cc)*ef(i,j+1) + &
951                (gg-dd)*ef(i+1,j+1) ) / (3.0*gg)
952       ENDIF
953
954    END SUBROUTINE collision_efficiency_rogers
955
[825]956 END MODULE lpm_collision_kernels_mod
Note: See TracBrowser for help on using the repository browser.