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

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

last commit documented

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