Line data Source code
1 : !===============================================================================
2 : ! Seasalt for Bulk Aerosol Model
3 : !===============================================================================
4 : module seasalt_model
5 : use shr_kind_mod, only: r8 => shr_kind_r8, cl => shr_kind_cl
6 : use ppgrid, only: pcols, pver
7 :
8 : implicit none
9 : private
10 :
11 : public :: seasalt_nbin
12 : public :: seasalt_nnum
13 : public :: seasalt_names
14 : public :: seasalt_indices
15 : public :: seasalt_init
16 : public :: seasalt_emis
17 : public :: seasalt_active
18 :
19 : public :: seasalt_depvel
20 :
21 : logical :: seasalt_active = .false.
22 :
23 : integer, parameter :: seasalt_nbin = 4
24 : integer, parameter :: seasalt_nnum = 0
25 :
26 : character(len=6), parameter :: seasalt_names(seasalt_nbin) &
27 : = (/'SSLT01', 'SSLT02', 'SSLT03', 'SSLT04'/)
28 :
29 : integer :: seasalt_indices(seasalt_nbin)
30 :
31 : contains
32 :
33 : !=============================================================================
34 : !=============================================================================
35 1536 : subroutine seasalt_init
36 : use cam_history, only: addfld, fieldname_len
37 : use constituents, only: cnst_get_ind
38 :
39 : character(len=fieldname_len) :: dummy
40 : integer :: m
41 :
42 7680 : do m = 1, seasalt_nbin
43 7680 : call cnst_get_ind(seasalt_names(m), seasalt_indices(m),abort=.false.)
44 : enddo
45 7680 : seasalt_active = any(seasalt_indices(:) > 0)
46 :
47 1536 : if (.not.seasalt_active) return
48 :
49 0 : dummy = 'RH'
50 0 : call addfld (dummy,(/ 'lev' /), 'A','frac','RH in dry dep calc')
51 0 : do m = 1,seasalt_nbin
52 0 : dummy = trim(seasalt_names(m)) // 'DI'
53 0 : call addfld (dummy,(/ 'lev' /), 'A','m/s',trim(seasalt_names(m))//' deposition diameter')
54 : enddo
55 :
56 1536 : end subroutine seasalt_init
57 :
58 : !=============================================================================
59 : !=============================================================================
60 0 : subroutine seasalt_emis( u10cubed, srf_temp, ocnfrc, ncol, cflx )
61 :
62 : ! dummy arguments
63 : real(r8), intent(in) :: u10cubed(:)
64 : real(r8), intent(in) :: srf_temp(:)
65 : real(r8), intent(in) :: ocnfrc(:)
66 : integer, intent(in) :: ncol
67 : real(r8), intent(inout) :: cflx(:,:)
68 :
69 : ! local vars
70 : integer :: ix,m
71 : real(r8), parameter :: sslt_source(seasalt_nbin) = (/ 4.77e-15_r8, 5.19e-14_r8, 1.22e-13_r8, 6.91e-14_r8 /)
72 :
73 0 : do m = 1, seasalt_nbin
74 0 : ix = seasalt_indices(m)
75 0 : cflx(:ncol,ix) = sslt_source(m) * u10cubed(:ncol) * ocnfrc(:ncol)
76 : enddo
77 :
78 1536 : end subroutine seasalt_emis
79 :
80 : !=============================================================================
81 : !=============================================================================
82 0 : subroutine seasalt_depvel( temp, pmid, q, ram1, fv, ncol, lchnk, vlc_dry,vlc_trb,vlc_grv )
83 : use aerosol_depvel, only: aerosol_depvel_compute
84 : use wv_saturation, only: qsat
85 : use cam_history, only: outfld
86 : use mo_constants, only: dns_aer_sst=>seasalt_density
87 :
88 : integer, intent(in) :: ncol, lchnk
89 : real(r8), intent(in) :: temp(:,:) ! temperature
90 : real(r8), intent(in) :: pmid(:,:) ! mid point pressure
91 : real(r8), intent(in) :: q(:,:) ! water vapor
92 : real(r8), intent(in) :: ram1(:) ! aerodynamical resistance (s/m)
93 : real(r8), intent(in) :: fv(:) ! friction velocity (m/s)
94 :
95 : real(r8), intent(out) :: vlc_trb(:,:) !Turbulent deposn velocity (m/s)
96 : real(r8), intent(out) :: vlc_grv(:,:,:) !grav deposn velocity (m/s)
97 : real(r8), intent(out) :: vlc_dry(:,:,:) !dry deposn velocity (m/s)
98 :
99 : real(r8) :: r
100 : real(r8) :: wetdia(pcols,pver,seasalt_nbin)
101 : real(r8) :: RH(pcols,pver),es(pcols,pver),qs(pcols,pver) ! for wet radius calculation
102 : real(r8),parameter:: c1=0.7674_r8, c2=3.0790_r8, c3=2.57e-11_r8,c4=-1.424_r8 ! wet radius calculation constants
103 : integer :: m, i,k
104 :
105 : ! set stokes correction to 1.0 for now not a bad assumption for our size range)
106 : real(r8), parameter :: sslt_stk_crc(seasalt_nbin) = (/ 1.0_r8, 1.0_r8, 1.0_r8, 1.0_r8 /)
107 : real(r8), parameter :: sslt_smt_vwr(seasalt_nbin) = (/0.52e-6_r8,2.38e-6_r8,4.86e-6_r8,15.14e-6_r8/)
108 :
109 : !-----------------------------------------------------------------------
110 0 : do k = 1, pver
111 0 : call qsat(temp(1:ncol,k),pmid(1:ncol,k),es(1:ncol,k),qs(1:ncol,k),ncol)
112 : end do
113 0 : RH(:ncol,:)=q(:ncol,:)/qs(:ncol,:)
114 0 : RH(:ncol,:)=max(0.01_r8,min(0.99_r8,RH(:ncol,:)))
115 : ! set stokes correction to 1.0 for now not a bad assumption for our size range)
116 0 : do m=1,seasalt_nbin
117 0 : r=sslt_smt_vwr(m)/2.0_r8
118 0 : do k=1,pver
119 0 : do i=1,ncol
120 0 : wetdia(i,k,m)=((r**3+c1*r**c2/(c3*r**c4-log(RH(i,k))))**(1._r8/3._r8))*2.0_r8
121 : enddo
122 : enddo
123 0 : call outfld( trim(seasalt_names(m))//'DI',wetdia(:,:,m), pcols, lchnk)
124 : enddo
125 0 : call outfld( 'RH',RH(:,:), pcols, lchnk)
126 :
127 : call aerosol_depvel_compute( ncol, pver, seasalt_nbin, temp, pmid, ram1, fv, wetdia, sslt_stk_crc, dns_aer_sst, &
128 0 : vlc_dry,vlc_trb,vlc_grv)
129 :
130 0 : endsubroutine seasalt_depvel
131 :
132 : end module seasalt_model
|