MODULE lpm_collision_kernels_mod !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: lpm_collision_kernels.f90 1093 2013-02-02 12:58:49Z raasch $ ! ! 1092 2013-02-02 11:24:22Z raasch ! unused variables removed ! ! 1071 2012-11-29 16:54:55Z franke ! Bugfix: collision efficiencies for Hall kernel should not be < 1.0E-20 ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1007 2012-09-19 14:30:36Z franke ! converted all units to SI units and replaced some parameters by corresponding ! PALM parameters ! Bugfix: factor in calculation of enhancement factor for collision efficencies ! changed from 10. to 1.0 ! ! 849 2012-03-15 10:35:09Z raasch ! routine collision_efficiency_rogers added (moved from former advec_particles ! to here) ! ! 835 2012-02-22 11:21:19Z raasch $ ! Bugfix: array diss can be used only in case of Wang kernel ! ! 828 2012-02-21 12:00:36Z raasch ! code has been completely reformatted, routine colker renamed ! recalculate_kernel, ! routine init_kernels added, radius is now communicated to the collision ! routines by array radclass ! ! Bugfix: transformation factor for dissipation changed from 1E5 to 1E4 ! ! 825 2012-02-19 03:03:44Z raasch ! routine renamed from wang_kernel to lpm_collision_kernels, ! turbulence_effects on collision replaced by wang_kernel ! ! 799 2011-12-21 17:48:03Z franke ! speed optimizations and formatting ! Bugfix: iq=1 is not allowed (routine effic) ! Bugfix: replaced stop by ec=0.0 in case of very small ec (routine effic) ! ! 790 2011-11-29 03:11:20Z raasch ! initial revision ! ! Description: ! ------------ ! This module calculates collision efficiencies either due to pure gravitational ! effects (Hall kernel, see Hall, 1980: J. Atmos. Sci., 2486-2507) or ! including the effects of (SGS) turbulence (Wang kernel, see Wang and ! Grabowski, 2009: Atmos. Sci. Lett., 10, 1-8). The original code has been ! provided by L.-P. Wang but is substantially reformatted and speed optimized ! here. ! ! ATTENTION: ! Physical quantities (like g, densities, etc.) used in this module still ! have to be adjusted to those values used in the main PALM code. ! Also, quantities in CGS-units should be converted to SI-units eventually. !------------------------------------------------------------------------------! USE arrays_3d USE cloud_parameters USE constants USE particle_attributes USE pegrid IMPLICIT NONE PRIVATE PUBLIC ckernel, collision_efficiency_rogers, init_kernels, & rclass_lbound, rclass_ubound, recalculate_kernel REAL :: epsilon, eps2, rclass_lbound, rclass_ubound, urms, urms2 REAL, DIMENSION(:), ALLOCATABLE :: epsclass, radclass, winf REAL, DIMENSION(:,:), ALLOCATABLE :: ec, ecf, gck, hkernel, hwratio REAL, DIMENSION(:,:,:), ALLOCATABLE :: ckernel SAVE ! !-- Public interfaces INTERFACE collision_efficiency_rogers MODULE PROCEDURE collision_efficiency_rogers END INTERFACE collision_efficiency_rogers INTERFACE init_kernels MODULE PROCEDURE init_kernels END INTERFACE init_kernels INTERFACE recalculate_kernel MODULE PROCEDURE recalculate_kernel END INTERFACE recalculate_kernel CONTAINS SUBROUTINE init_kernels !------------------------------------------------------------------------------! ! Initialization of the collision efficiency matrix with fixed radius and ! dissipation classes, calculated at simulation start only. !------------------------------------------------------------------------------! IMPLICIT NONE INTEGER :: i, j, k ! !-- Calculate collision efficiencies for fixed radius- and dissipation !-- classes IF ( collision_kernel(6:9) == 'fast' ) THEN ALLOCATE( ckernel(1:radius_classes,1:radius_classes, & 0:dissipation_classes), epsclass(1:dissipation_classes), & radclass(1:radius_classes) ) ! !-- Calculate the radius class bounds with logarithmic distances !-- in the interval [1.0E-6, 2.0E-4] m rclass_lbound = LOG( 1.0E-6 ) rclass_ubound = LOG( 2.0E-4 ) radclass(1) = 1.0E-6 DO i = 2, radius_classes radclass(i) = EXP( rclass_lbound + & ( rclass_ubound - rclass_lbound ) * ( i-1.0 ) /& ( radius_classes - 1.0 ) ) ! IF ( myid == 0 ) THEN ! PRINT*, 'i=', i, ' r = ', radclass(i)*1.0E6 ! ENDIF ENDDO ! !-- Set the class bounds for dissipation in interval [0.0, 0.1] m**2/s**3 DO i = 1, dissipation_classes epsclass(i) = 0.1 * REAL( i ) / dissipation_classes ! IF ( myid == 0 ) THEN ! PRINT*, 'i=', i, ' eps = ', epsclass(i) ! ENDIF ENDDO ! !-- Calculate collision efficiencies of the Wang/ayala kernel ALLOCATE( ec(1:radius_classes,1:radius_classes), & ecf(1:radius_classes,1:radius_classes), & gck(1:radius_classes,1:radius_classes), & winf(1:radius_classes) ) DO k = 1, dissipation_classes epsilon = epsclass(k) urms = 2.02 * ( epsilon / 0.04 )**( 1.0 / 3.0 ) CALL turbsd CALL turb_enhance_eff CALL effic DO j = 1, radius_classes DO i = 1, radius_classes ckernel(i,j,k) = ec(i,j) * gck(i,j) * ecf(i,j) ENDDO ENDDO ENDDO ! !-- Calculate collision efficiencies of the Hall kernel ALLOCATE( hkernel(1:radius_classes,1:radius_classes), & hwratio(1:radius_classes,1:radius_classes) ) CALL fallg CALL effic DO j = 1, radius_classes DO i = 1, radius_classes hkernel(i,j) = pi * ( radclass(j) + radclass(i) )**2 & * ec(i,j) * ABS( winf(j) - winf(i) ) ckernel(i,j,0) = hkernel(i,j) ! hall kernel stored on index 0 ENDDO ENDDO ! !-- Test output of efficiencies IF ( j == -1 ) THEN PRINT*, '*** Hall kernel' WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i)*1.0E6, & i = 1,radius_classes ) DO j = 1, radius_classes WRITE ( *,'(F4.0,1X,20(F8.4,1X))' ) radclass(j), & ( hkernel(i,j), i = 1,radius_classes ) ENDDO DO k = 1, dissipation_classes DO i = 1, radius_classes DO j = 1, radius_classes IF ( hkernel(i,j) == 0.0 ) THEN hwratio(i,j) = 9999999.9 ELSE hwratio(i,j) = ckernel(i,j,k) / hkernel(i,j) ENDIF ENDDO ENDDO PRINT*, '*** epsilon = ', epsclass(k) WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i)*1.0E6, & i = 1,radius_classes ) DO j = 1, radius_classes ! WRITE ( *,'(F4.0,1X,20(F4.2,1X))' ) radclass(j)*1.0E6, & ! ( ckernel(i,j,k), i = 1,radius_classes ) WRITE ( *,'(F4.0,1X,20(F8.4,1X))' ) radclass(j)*1.0E6, & ( hwratio(i,j), i = 1,radius_classes ) ENDDO ENDDO ENDIF DEALLOCATE( ec, ecf, epsclass, gck, hkernel, winf ) ELSEIF( collision_kernel == 'hall' .OR. collision_kernel == 'wang' ) & THEN ! !-- Initial settings for Hall- and Wang-Kernel !-- To be done: move here parts from turbsd, fallg, ecoll, etc. ENDIF END SUBROUTINE init_kernels !------------------------------------------------------------------------------! ! Calculation of collision kernels during each timestep and for each grid box !------------------------------------------------------------------------------! SUBROUTINE recalculate_kernel( i1, j1, k1 ) USE arrays_3d USE cloud_parameters USE constants USE cpulog USE indices USE interfaces USE particle_attributes IMPLICIT NONE INTEGER :: i, i1, j, j1, k1, pend, pstart pstart = prt_start_index(k1,j1,i1) pend = prt_start_index(k1,j1,i1) + prt_count(k1,j1,i1) - 1 radius_classes = prt_count(k1,j1,i1) ALLOCATE( ec(1:radius_classes,1:radius_classes), & radclass(1:radius_classes), winf(1:radius_classes) ) ! !-- Store particle radii on the radclass array radclass(1:radius_classes) = particles(pstart:pend)%radius IF ( wang_kernel ) THEN epsilon = diss(k1,j1,i1) ! dissipation rate in m**2/s**3 ELSE epsilon = 0.0 ENDIF urms = 2.02 * ( epsilon / 0.04 )**( 0.33333333333 ) IF ( wang_kernel .AND. epsilon > 1.0E-7 ) THEN ! !-- Call routines to calculate efficiencies for the Wang kernel ALLOCATE( gck(1:radius_classes,1:radius_classes), & ecf(1:radius_classes,1:radius_classes) ) CALL turbsd CALL turb_enhance_eff CALL effic DO j = 1, radius_classes DO i = 1, radius_classes ckernel(pstart+i-1,pstart+j-1,1) = ec(i,j) * gck(i,j) * ecf(i,j) ENDDO ENDDO DEALLOCATE( gck, ecf ) ELSE ! !-- Call routines to calculate efficiencies for the Hall kernel CALL fallg CALL effic DO j = 1, radius_classes DO i = 1, radius_classes ckernel(pstart+i-1,pstart+j-1,1) = pi * & ( radclass(j) + radclass(i) )**2 & * ec(i,j) * ABS( winf(j) - winf(i) ) ENDDO ENDDO ENDIF DEALLOCATE( ec, radclass, winf ) END SUBROUTINE recalculate_kernel !------------------------------------------------------------------------------! ! Calculation of gck ! This is from Aayala 2008b, page 37ff. ! Necessary input parameters: water density, radii of droplets, air density, ! air viscosity, turbulent dissipation rate, taylor microscale reynolds number, ! gravitational acceleration --> to be replaced by PALM parameters !------------------------------------------------------------------------------! SUBROUTINE turbsd USE constants USE cloud_parameters USE particle_attributes USE arrays_3d USE control_parameters IMPLICIT NONE INTEGER :: i, j LOGICAL, SAVE :: first = .TRUE. REAL :: ao, ao_gr, bbb, be, b1, b2, ccc, c1, c1_gr, c2, d1, d2, eta, & e1, e2, fao_gr, fr, grfin, lambda, lambda_re, lf, rc, rrp, & sst, tauk, tl, t2, tt, t1, vk, vrms1xy, vrms2xy, v1, v1v2xy, & v1xysq, v2, v2xysq, wrfin, wrgrav2, wrtur2xy, xx, yy, z REAL, DIMENSION(1:radius_classes) :: st, tau ! !-- Initial assignment of constants IF ( first ) THEN first = .FALSE. ENDIF lambda = urms * SQRT( 15.0 * molecular_viscosity / epsilon ) ! in m lambda_re = urms**2 * SQRT( 15.0 / epsilon / molecular_viscosity ) tl = urms**2 / epsilon ! in s lf = 0.5 * urms**3 / epsilon ! in m tauk = SQRT( molecular_viscosity / epsilon ) ! in s eta = ( molecular_viscosity**3 / epsilon )**0.25 ! in m vk = eta / tauk ao = ( 11.0 + 7.0 * lambda_re ) / ( 205.0 + lambda_re ) tt = SQRT( 2.0 * lambda_re / ( SQRT( 15.0 ) * ao ) ) * tauk ! in s CALL fallg ! gives winf in m/s DO i = 1, radius_classes tau(i) = winf(i) / g ! in s st(i) = tau(i) / tauk ENDDO ! !-- Calculate wr (from Aayala 2008b, page 38f) z = tt / tl be = SQRT( 2.0 ) * lambda / lf bbb = SQRT( 1.0 - 2.0 * be**2 ) d1 = ( 1.0 + bbb ) / ( 2.0 * bbb ) e1 = lf * ( 1.0 + bbb ) * 0.5 ! in m d2 = ( 1.0 - bbb ) * 0.5 / bbb e2 = lf * ( 1.0 - bbb ) * 0.5 ! in m ccc = SQRT( 1.0 - 2.0 * z**2 ) b1 = ( 1.0 + ccc ) * 0.5 / ccc c1 = tl * ( 1.0 + ccc ) * 0.5 ! in s b2 = ( 1.0 - ccc ) * 0.5 / ccc c2 = tl * ( 1.0 - ccc ) * 0.5 ! in s DO i = 1, radius_classes v1 = winf(i) ! in m/s t1 = tau(i) ! in s DO j = 1, i rrp = radclass(i) + radclass(j) v2 = winf(j) ! in m/s t2 = tau(j) ! in s v1xysq = b1 * d1 * phi_w(c1,e1,v1,t1) - b1 * d2 * phi_w(c1,e2,v1,t1) & - b2 * d1 * phi_w(c2,e1,v1,t1) + b2 * d2 * phi_w(c2,e2,v1,t1) v1xysq = v1xysq * urms**2 / t1 ! in m**2/s**2 vrms1xy = SQRT( v1xysq ) ! in m/s v2xysq = b1 * d1 * phi_w(c1,e1,v2,t2) - b1 * d2 * phi_w(c1,e2,v2,t2) & - b2 * d1 * phi_w(c2,e1,v2,t2) + b2 * d2 * phi_w(c2,e2,v2,t2) v2xysq = v2xysq * urms**2 / t2 ! in m**2/s**2 vrms2xy = SQRT( v2xysq ) ! in m/s IF ( winf(i) >= winf(j) ) THEN v1 = winf(i) t1 = tau(i) v2 = winf(j) t2 = tau(j) ELSE v1 = winf(j) t1 = tau(j) v2 = winf(i) t2 = tau(i) ENDIF v1v2xy = b1 * d1 * zhi(c1,e1,v1,t1,v2,t2) - & b1 * d2 * zhi(c1,e2,v1,t1,v2,t2) - & b2 * d1 * zhi(c2,e1,v1,t1,v2,t2) + & b2 * d2* zhi(c2,e2,v1,t1,v2,t2) fr = d1 * EXP( -rrp / e1 ) - d2 * EXP( -rrp / e2 ) v1v2xy = v1v2xy * fr * urms**2 / tau(i) / tau(j) ! in m**2/s**2 wrtur2xy = vrms1xy**2 + vrms2xy**2 - 2.0 * v1v2xy ! in m**2/s**2 IF ( wrtur2xy < 0.0 ) wrtur2xy = 0.0 wrgrav2 = pi / 8.0 * ( winf(j) - winf(i) )**2 wrfin = SQRT( ( 2.0 / pi ) * ( wrtur2xy + wrgrav2) ) ! in m/s ! !-- Calculate gr IF ( st(j) > st(i) ) THEN sst = st(j) ELSE sst = st(i) ENDIF xx = -0.1988 * sst**4 + 1.5275 * sst**3 - 4.2942 * sst**2 + & 5.3406 * sst IF ( xx < 0.0 ) xx = 0.0 yy = 0.1886 * EXP( 20.306 / lambda_re ) c1_gr = xx / ( g / vk * tauk )**yy ao_gr = ao + ( pi / 8.0) * ( g / vk * tauk )**2 fao_gr = 20.115 * SQRT( ao_gr / lambda_re ) rc = SQRT( fao_gr * ABS( st(j) - st(i) ) ) * eta ! in cm grfin = ( ( eta**2 + rc**2 ) / ( rrp**2 + rc**2) )**( c1_gr*0.5 ) IF ( grfin < 1.0 ) grfin = 1.0 gck(i,j) = 2.0 * pi * rrp**2 * wrfin * grfin ! in cm**3/s gck(j,i) = gck(i,j) ENDDO ENDDO END SUBROUTINE turbsd !------------------------------------------------------------------------------! ! phi_w as a function !------------------------------------------------------------------------------! REAL FUNCTION phi_w( a, b, vsett, tau0 ) IMPLICIT NONE REAL :: a, aa1, b, tau0, vsett aa1 = 1.0 / tau0 + 1.0 / a + vsett / b phi_w = 1.0 / aa1 - 0.5 * vsett / b / aa1**2 ! in s END FUNCTION phi_w !------------------------------------------------------------------------------! ! zhi as a function !------------------------------------------------------------------------------! REAL FUNCTION zhi( a, b, vsett1, tau1, vsett2, tau2 ) IMPLICIT NONE REAL :: a, aa1, aa2, aa3, aa4, aa5, aa6, b, tau1, tau2, vsett1, vsett2 aa1 = vsett2 / b - 1.0 / tau2 - 1.0 / a aa2 = vsett1 / b + 1.0 / tau1 + 1.0 / a aa3 = ( vsett1 - vsett2 ) / b + 1.0 / tau1 + 1.0 / tau2 aa4 = ( vsett2 / b )**2 - ( 1.0 / tau2 + 1.0 / a )**2 aa5 = vsett2 / b + 1.0 / tau2 + 1.0 / a aa6 = 1.0 / tau1 - 1.0 / a + ( 1.0 / tau2 + 1.0 / a) * vsett1 / vsett2 zhi = (1.0 / aa1 - 1.0 / aa2 ) * ( vsett1 - vsett2 ) * 0.5 / b / aa3**2 & + (4.0 / aa4 - 1.0 / aa5**2 - 1.0 / aa1**2 ) * vsett2 * 0.5 / b /aa6& + (2.0 * ( b / aa2 - b / aa1 ) - vsett1 / aa2**2 + vsett2 / aa1**2 )& * 0.5 / b / aa3 ! in s**2 END FUNCTION zhi !------------------------------------------------------------------------------! ! Calculation of terminal velocity winf following Equations 10-138 to 10-145 ! from (Pruppacher and Klett, 1997) !------------------------------------------------------------------------------! SUBROUTINE fallg USE constants USE cloud_parameters USE particle_attributes USE arrays_3d USE control_parameters IMPLICIT NONE INTEGER :: i, j LOGICAL, SAVE :: first = .TRUE. REAL, SAVE :: cunh, eta, phy, py, rho_a, sigma, stb, stok, xlamb REAL :: bond, x, xrey, y REAL, DIMENSION(1:7), SAVE :: b REAL, DIMENSION(1:6), SAVE :: c ! !-- Initial assignment of constants IF ( first ) THEN first = .FALSE. b = (/ -0.318657E1, 0.992696E0, -0.153193E-2, -0.987059E-3, & -0.578878E-3, 0.855176E-4, -0.327815E-5 /) c = (/ -0.500015E1, 0.523778E1, -0.204914E1, 0.475294E0, & -0.542819E-1, 0.238449E-2 /) ! !-- Parameter values for p = 1013,25 hPa and T = 293,15 K eta = 1.818E-5 ! in kg/(m s) xlamb = 6.6E-8 ! in m rho_a = 1.204 ! in kg/m**3 cunh = 1.26 * xlamb ! in m sigma = 0.07363 ! in kg/s**2 stok = 2.0 * g * ( rho_l - rho_a ) / ( 9.0 * eta ) ! in 1/(m s) stb = 32.0 * rho_a * ( rho_l - rho_a) * g / (3.0 * eta * eta) phy = sigma**3 * rho_a**2 / ( eta**4 * g * ( rho_l - rho_a ) ) py = phy**( 1.0 / 6.0 ) ENDIF DO j = 1, radius_classes IF ( radclass(j) <= 1.0E-5 ) THEN winf(j) = stok * ( radclass(j)**2 + cunh * radclass(j) ) ELSEIF ( radclass(j) > 1.0E-5 .AND. radclass(j) <= 5.35E-4 ) THEN x = LOG( stb * radclass(j)**3 ) y = 0.0 DO i = 1, 7 y = y + b(i) * x**(i-1) ENDDO ! !-- Note: this Eq. is wrong in (Pruppacher and Klett, 1997, p. 418) !-- for correct version see (Beard, 1976) xrey = ( 1.0 + cunh / radclass(j) ) * EXP( y ) winf(j) = xrey * eta / ( 2.0 * rho_a * radclass(j) ) ELSEIF ( radclass(j) > 5.35E-4 ) THEN IF ( radclass(j) > 0.0035 ) THEN bond = g * ( rho_l - rho_a ) * 0.0035**2 / sigma ELSE bond = g * ( rho_l - rho_a ) * radclass(j)**2 / sigma ENDIF x = LOG( 16.0 * bond * py / 3.0 ) y = 0.0 DO i = 1, 6 y = y + c(i) * x**(i-1) ENDDO xrey = py * EXP( y ) IF ( radclass(j) > 0.0035 ) THEN winf(j) = xrey * eta / ( 2.0 * rho_a * 0.0035 ) ELSE winf(j) = xrey * eta / ( 2.0 * rho_a * radclass(j) ) ENDIF ENDIF ENDDO END SUBROUTINE fallg !------------------------------------------------------------------------------! ! Calculation of collision efficiencies for the Hall kernel !------------------------------------------------------------------------------! SUBROUTINE effic USE arrays_3d USE cloud_parameters USE constants USE particle_attributes IMPLICIT NONE INTEGER :: i, iq, ir, j, k INTEGER, DIMENSION(:), ALLOCATABLE :: ira LOGICAL, SAVE :: first = .TRUE. REAL :: ek, particle_radius, pp, qq, rq REAL, DIMENSION(1:21), SAVE :: rat REAL, DIMENSION(1:15), SAVE :: r0 REAL, DIMENSION(1:15,1:21), SAVE :: ecoll ! !-- Initial assignment of constants IF ( first ) THEN first = .FALSE. r0 = (/ 6.0, 8.0, 10.0, 15.0, 20.0, 25.0, 30.0, 40.0, 50.0, 60., & 70.0, 100.0, 150.0, 200.0, 300.0 /) rat = (/ 0.00, 0.05, 0.10, 0.15, 0.20, 0.25, 0.30, 0.35, 0.40, 0.45, & 0.50, 0.55, 0.60, 0.65, 0.70, 0.75, 0.80, 0.85, 0.90, 0.95, & 1.00 /) ecoll(:,1) = (/0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, & 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001/) ecoll(:,2) = (/0.003, 0.003, 0.003, 0.004, 0.005, 0.005, 0.005, & 0.010, 0.100, 0.050, 0.200, 0.500, 0.770, 0.870, 0.970/) ecoll(:,3) = (/0.007, 0.007, 0.007, 0.008, 0.009, 0.010, 0.010, & 0.070, 0.400, 0.430, 0.580, 0.790, 0.930, 0.960, 1.000/) ecoll(:,4) = (/0.009, 0.009, 0.009, 0.012, 0.015, 0.010, 0.020, & 0.280, 0.600, 0.640, 0.750, 0.910, 0.970, 0.980, 1.000/) ecoll(:,5) = (/0.014, 0.014, 0.014, 0.015, 0.016, 0.030, 0.060, & 0.500, 0.700, 0.770, 0.840, 0.950, 0.970, 1.000, 1.000/) ecoll(:,6) = (/0.017, 0.017, 0.017, 0.020, 0.022, 0.060, 0.100, & 0.620, 0.780, 0.840, 0.880, 0.950, 1.000, 1.000, 1.000/) ecoll(:,7) = (/0.030, 0.030, 0.024, 0.022, 0.032, 0.062, 0.200, & 0.680, 0.830, 0.870, 0.900, 0.950, 1.000, 1.000, 1.000/) ecoll(:,8) = (/0.025, 0.025, 0.025, 0.036, 0.043, 0.130, 0.270, & 0.740, 0.860, 0.890, 0.920, 1.000, 1.000, 1.000, 1.000/) ecoll(:,9) = (/0.027, 0.027, 0.027, 0.040, 0.052, 0.200, 0.400, & 0.780, 0.880, 0.900, 0.940, 1.000, 1.000, 1.000, 1.000/) ecoll(:,10)= (/0.030, 0.030, 0.030, 0.047, 0.064, 0.250, 0.500, & 0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,11)= (/0.040, 0.040, 0.033, 0.037, 0.068, 0.240, 0.550, & 0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,12)= (/0.035, 0.035, 0.035, 0.055, 0.079, 0.290, 0.580, & 0.800, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,13)= (/0.037, 0.037, 0.037, 0.062, 0.082, 0.290, 0.590, & 0.780, 0.900, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,14)= (/0.037, 0.037, 0.037, 0.060, 0.080, 0.290, 0.580, & 0.770, 0.890, 0.910, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,15)= (/0.037, 0.037, 0.037, 0.041, 0.075, 0.250, 0.540, & 0.760, 0.880, 0.920, 0.950, 1.000, 1.000, 1.000, 1.000/) ecoll(:,16)= (/0.037, 0.037, 0.037, 0.052, 0.067, 0.250, 0.510, & 0.770, 0.880, 0.930, 0.970, 1.000, 1.000, 1.000, 1.000/) ecoll(:,17)= (/0.037, 0.037, 0.037, 0.047, 0.057, 0.250, 0.490, & 0.770, 0.890, 0.950, 1.000, 1.000, 1.000, 1.000, 1.000/) ecoll(:,18)= (/0.036, 0.036, 0.036, 0.042, 0.048, 0.230, 0.470, & 0.780, 0.920, 1.000, 1.020, 1.020, 1.020, 1.020, 1.020/) ecoll(:,19)= (/0.040, 0.040, 0.035, 0.033, 0.040, 0.112, 0.450, & 0.790, 1.010, 1.030, 1.040, 1.040, 1.040, 1.040, 1.040/) ecoll(:,20)= (/0.033, 0.033, 0.033, 0.033, 0.033, 0.119, 0.470, & 0.950, 1.300, 1.700, 2.300, 2.300, 2.300, 2.300, 2.300/) ecoll(:,21)= (/0.027, 0.027, 0.027, 0.027, 0.027, 0.125, 0.520, & 1.400, 2.300, 3.000, 4.000, 4.000, 4.000, 4.000, 4.000/) ENDIF ! !-- Calculate the radius class index of particles with respect to array r !-- Radius has to be in µm ALLOCATE( ira(1:radius_classes) ) DO j = 1, radius_classes particle_radius = radclass(j) * 1.0E6 DO k = 1, 15 IF ( particle_radius < r0(k) ) THEN ira(j) = k EXIT ENDIF ENDDO IF ( particle_radius >= r0(15) ) ira(j) = 16 ENDDO ! !-- Two-dimensional linear interpolation of the collision efficiency. !-- Radius has to be in µm DO j = 1, radius_classes DO i = 1, j ir = ira(j) rq = radclass(i) / radclass(j) iq = INT( rq * 20 ) + 1 iq = MAX( iq , 2) IF ( ir < 16 ) THEN IF ( ir >= 2 ) THEN pp = ( ( radclass(j) * 1.0E06 ) - r0(ir-1) ) / & ( r0(ir) - r0(ir-1) ) qq = ( rq- rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) ec(j,i) = ( 1.0-pp ) * ( 1.0-qq ) * ecoll(ir-1,iq-1) & + pp * ( 1.0-qq ) * ecoll(ir,iq-1) & + qq * ( 1.0-pp ) * ecoll(ir-1,iq) & + pp * qq * ecoll(ir,iq) ELSE qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) ec(j,i) = (1.0-qq) * ecoll(1,iq-1) + qq * ecoll(1,iq) ENDIF ELSE qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) ek = ( 1.0 - qq ) * ecoll(15,iq-1) + qq * ecoll(15,iq) ec(j,i) = MIN( ek, 1.0 ) ENDIF IF ( ec(j,i) < 1.0E-20 ) ec(j,i) = 0.0 ec(i,j) = ec(j,i) ENDDO ENDDO DEALLOCATE( ira ) END SUBROUTINE effic !------------------------------------------------------------------------------! ! Calculation of enhancement factor for collision efficencies due to turbulence !------------------------------------------------------------------------------! SUBROUTINE turb_enhance_eff USE constants USE cloud_parameters USE particle_attributes USE arrays_3d IMPLICIT NONE INTEGER :: i, iq, ir, j, k, kk INTEGER, DIMENSION(:), ALLOCATABLE :: ira REAL :: particle_radius, pp, qq, rq, y1, y2, y3 LOGICAL, SAVE :: first = .TRUE. REAL, DIMENSION(1:11), SAVE :: rat REAL, DIMENSION(1:7), SAVE :: r0 REAL, DIMENSION(1:7,1:11), SAVE :: ecoll_100, ecoll_400 ! !-- Initial assignment of constants IF ( first ) THEN first = .FALSE. r0 = (/ 10.0, 20.0, 30.0, 40.0, 50.0, 60.0, 100.0 /) rat = (/ 0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0 /) ! !-- for 100 cm**2/s**3 ecoll_100(:,1) = (/1.74, 1.74, 1.773, 1.49, 1.207, 1.207, 1.0 /) ecoll_100(:,2) = (/1.46, 1.46, 1.421, 1.245, 1.069, 1.069, 1.0 /) ecoll_100(:,3) = (/1.32, 1.32, 1.245, 1.123, 1.000, 1.000, 1.0 /) ecoll_100(:,4) = (/1.250, 1.250, 1.148, 1.087, 1.025, 1.025, 1.0 /) ecoll_100(:,5) = (/1.186, 1.186, 1.066, 1.060, 1.056, 1.056, 1.0 /) ecoll_100(:,6) = (/1.045, 1.045, 1.000, 1.014, 1.028, 1.028, 1.0 /) ecoll_100(:,7) = (/1.070, 1.070, 1.030, 1.038, 1.046, 1.046, 1.0 /) ecoll_100(:,8) = (/1.000, 1.000, 1.054, 1.042, 1.029, 1.029, 1.0 /) ecoll_100(:,9) = (/1.223, 1.223, 1.117, 1.069, 1.021, 1.021, 1.0 /) ecoll_100(:,10)= (/1.570, 1.570, 1.244, 1.166, 1.088, 1.088, 1.0 /) ecoll_100(:,11)= (/20.3, 20.3, 14.6 , 8.61, 2.60, 2.60 , 1.0 /) ! !-- for 400 cm**2/s**3 ecoll_400(:,1) = (/4.976, 4.976, 3.593, 2.519, 1.445, 1.445, 1.0 /) ecoll_400(:,2) = (/2.984, 2.984, 2.181, 1.691, 1.201, 1.201, 1.0 /) ecoll_400(:,3) = (/1.988, 1.988, 1.475, 1.313, 1.150, 1.150, 1.0 /) ecoll_400(:,4) = (/1.490, 1.490, 1.187, 1.156, 1.126, 1.126, 1.0 /) ecoll_400(:,5) = (/1.249, 1.249, 1.088, 1.090, 1.092, 1.092, 1.0 /) ecoll_400(:,6) = (/1.139, 1.139, 1.130, 1.091, 1.051, 1.051, 1.0 /) ecoll_400(:,7) = (/1.220, 1.220, 1.190, 1.138, 1.086, 1.086, 1.0 /) ecoll_400(:,8) = (/1.325, 1.325, 1.267, 1.165, 1.063, 1.063, 1.0 /) ecoll_400(:,9) = (/1.716, 1.716, 1.345, 1.223, 1.100, 1.100, 1.0 /) ecoll_400(:,10)= (/3.788, 3.788, 1.501, 1.311, 1.120, 1.120, 1.0 /) ecoll_400(:,11)= (/36.52, 36.52, 19.16, 22.80, 26.0, 26.0, 1.0 /) ENDIF ! !-- Calculate the radius class index of particles with respect to array r0 !-- Radius has to be in µm ALLOCATE( ira(1:radius_classes) ) DO j = 1, radius_classes particle_radius = radclass(j) * 1.0E6 DO k = 1, 7 IF ( particle_radius < r0(k) ) THEN ira(j) = k EXIT ENDIF ENDDO IF ( particle_radius >= r0(7) ) ira(j) = 8 ENDDO ! !-- Two-dimensional linear interpolation of the collision efficiencies !-- Radius has to be in µm DO j = 1, radius_classes DO i = 1, j ir = ira(j) rq = radclass(i) / radclass(j) DO kk = 2, 11 IF ( rq <= rat(kk) ) THEN iq = kk EXIT ENDIF ENDDO y1 = 0.0001 ! for 0 m**2/s**3 IF ( ir < 8 ) THEN IF ( ir >= 2 ) THEN pp = ( radclass(j)*1.0E6 - r0(ir-1) ) / ( r0(ir) - r0(ir-1) ) qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) y2 = ( 1.0-pp ) * ( 1.0-qq ) * ecoll_100(ir-1,iq-1) + & pp * ( 1.0-qq ) * ecoll_100(ir,iq-1) + & qq * ( 1.0-pp ) * ecoll_100(ir-1,iq) + & pp * qq * ecoll_100(ir,iq) y3 = ( 1.0-pp ) * ( 1.0-qq ) * ecoll_400(ir-1,iq-1) + & pp * ( 1.0-qq ) * ecoll_400(ir,iq-1) + & qq * ( 1.0-pp ) * ecoll_400(ir-1,iq) + & pp * qq * ecoll_400(ir,iq) ELSE qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) y2 = ( 1.0-qq ) * ecoll_100(1,iq-1) + qq * ecoll_100(1,iq) y3 = ( 1.0-qq ) * ecoll_400(1,iq-1) + qq * ecoll_400(1,iq) ENDIF ELSE qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) ) y2 = ( 1.0-qq ) * ecoll_100(7,iq-1) + qq * ecoll_100(7,iq) y3 = ( 1.0-qq ) * ecoll_400(7,iq-1) + qq * ecoll_400(7,iq) ENDIF ! !-- Linear interpolation of dissipation rate in m**2/s**3 IF ( epsilon <= 0.01 ) THEN ecf(j,i) = ( epsilon - 0.01 ) / ( 0.0 - 0.01 ) * y1 & + ( epsilon - 0.0 ) / ( 0.01 - 0.0 ) * y2 ELSEIF ( epsilon <= 0.06 ) THEN ecf(j,i) = ( epsilon - 0.04 ) / ( 0.01 - 0.04 ) * y2 & + ( epsilon - 0.01 ) / ( 0.04 - 0.01 ) * y3 ELSE ecf(j,i) = ( 0.06 - 0.04 ) / ( 0.01 - 0.04 ) * y2 & + ( 0.06 - 0.01 ) / ( 0.04 - 0.01 ) * y3 ENDIF IF ( ecf(j,i) < 1.0 ) ecf(j,i) = 1.0 ecf(i,j) = ecf(j,i) ENDDO ENDDO END SUBROUTINE turb_enhance_eff SUBROUTINE collision_efficiency_rogers( mean_r, r, e) !------------------------------------------------------------------------------! ! Collision efficiencies from table 8.2 in Rogers and Yau (1989, 3rd edition). ! Values are calculated from table by bilinear interpolation. !------------------------------------------------------------------------------! IMPLICIT NONE INTEGER :: i, j, k LOGICAL, SAVE :: first = .TRUE. REAL :: aa, bb, cc, dd, dx, dy, e, gg, mean_r, mean_rm, r, & rm, x, y REAL, DIMENSION(1:9), SAVE :: collected_r = 0.0 REAL, DIMENSION(1:19), SAVE :: collector_r = 0.0 REAL, DIMENSION(1:9,1:19), SAVE :: ef = 0.0 mean_rm = mean_r * 1.0E06 rm = r * 1.0E06 IF ( first ) THEN collected_r = (/ 2.0, 3.0, 4.0, 6.0, 8.0, 10.0, 15.0, 20.0, 25.0 /) collector_r = (/ 10.0, 20.0, 30.0, 40.0, 50.0, 60.0, 80.0, 100.0, & 150.0, 200.0, 300.0, 400.0, 500.0, 600.0, 1000.0, & 1400.0, 1800.0, 2400.0, 3000.0 /) ef(:,1) = (/0.017, 0.027, 0.037, 0.052, 0.052, 0.052, 0.052, 0.0, & 0.0 /) ef(:,2) = (/0.001, 0.016, 0.027, 0.060, 0.12, 0.17, 0.17, 0.17, 0.0 /) ef(:,3) = (/0.001, 0.001, 0.02, 0.13, 0.28, 0.37, 0.54, 0.55, 0.47/) ef(:,4) = (/0.001, 0.001, 0.02, 0.23, 0.4, 0.55, 0.7, 0.75, 0.75/) ef(:,5) = (/0.01, 0.01, 0.03, 0.3, 0.4, 0.58, 0.73, 0.75, 0.79/) ef(:,6) = (/0.01, 0.01, 0.13, 0.38, 0.57, 0.68, 0.80, 0.86, 0.91/) ef(:,7) = (/0.01, 0.085, 0.23, 0.52, 0.68, 0.76, 0.86, 0.92, 0.95/) ef(:,8) = (/0.01, 0.14, 0.32, 0.60, 0.73, 0.81, 0.90, 0.94, 0.96/) ef(:,9) = (/0.025, 0.25, 0.43, 0.66, 0.78, 0.83, 0.92, 0.95, 0.96/) ef(:,10)= (/0.039, 0.3, 0.46, 0.69, 0.81, 0.87, 0.93, 0.95, 0.96/) ef(:,11)= (/0.095, 0.33, 0.51, 0.72, 0.82, 0.87, 0.93, 0.96, 0.97/) ef(:,12)= (/0.098, 0.36, 0.51, 0.73, 0.83, 0.88, 0.93, 0.96, 0.97/) ef(:,13)= (/0.1, 0.36, 0.52, 0.74, 0.83, 0.88, 0.93, 0.96, 0.97/) ef(:,14)= (/0.17, 0.4, 0.54, 0.72, 0.83, 0.88, 0.94, 0.98, 1.0 /) ef(:,15)= (/0.15, 0.37, 0.52, 0.74, 0.82, 0.88, 0.94, 0.98, 1.0 /) ef(:,16)= (/0.11, 0.34, 0.49, 0.71, 0.83, 0.88, 0.94, 0.95, 1.0 /) ef(:,17)= (/0.08, 0.29, 0.45, 0.68, 0.8, 0.86, 0.96, 0.94, 1.0 /) ef(:,18)= (/0.04, 0.22, 0.39, 0.62, 0.75, 0.83, 0.92, 0.96, 1.0 /) ef(:,19)= (/0.02, 0.16, 0.33, 0.55, 0.71, 0.81, 0.90, 0.94, 1.0 /) ENDIF DO k = 1, 8 IF ( collected_r(k) <= mean_rm ) i = k ENDDO DO k = 1, 18 IF ( collector_r(k) <= rm ) j = k ENDDO IF ( rm < 10.0 ) THEN e = 0.0 ELSEIF ( mean_rm < 2.0 ) THEN e = 0.001 ELSEIF ( mean_rm >= 25.0 ) THEN IF( j <= 2 ) e = 0.0 IF( j == 3 ) e = 0.47 IF( j == 4 ) e = 0.8 IF( j == 5 ) e = 0.9 IF( j >=6 ) e = 1.0 ELSEIF ( rm >= 3000.0 ) THEN IF( i == 1 ) e = 0.02 IF( i == 2 ) e = 0.16 IF( i == 3 ) e = 0.33 IF( i == 4 ) e = 0.55 IF( i == 5 ) e = 0.71 IF( i == 6 ) e = 0.81 IF( i == 7 ) e = 0.90 IF( i >= 8 ) e = 0.94 ELSE x = mean_rm - collected_r(i) y = rm - collector_r(j) dx = collected_r(i+1) - collected_r(i) dy = collector_r(j+1) - collector_r(j) aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd e = ( (gg-aa)*ef(i,j) + (gg-bb)*ef(i+1,j) + (gg-cc)*ef(i,j+1) + & (gg-dd)*ef(i+1,j+1) ) / (3.0*gg) ENDIF END SUBROUTINE collision_efficiency_rogers END MODULE lpm_collision_kernels_mod