Line data Source code
1 : !===============================================================================
2 : ! Seasalt for Modal 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 : use modal_aero_data,only: ntot_amode, nslt=>nSeaSalt
8 :
9 : implicit none
10 : private
11 :
12 : public :: seasalt_nbin
13 : public :: seasalt_nnum
14 : public :: seasalt_names
15 : public :: seasalt_indices
16 : public :: seasalt_init
17 : public :: seasalt_emis
18 : public :: seasalt_active
19 :
20 : integer, protected :: seasalt_nbin ! = nslt
21 : integer, protected :: seasalt_nnum ! = nnum
22 :
23 : character(len=6), protected, allocatable :: seasalt_names(:)
24 : integer, protected, allocatable :: seasalt_indices(:)
25 :
26 : logical :: seasalt_active = .false.
27 :
28 : real(r8):: emis_scale
29 :
30 : contains
31 :
32 : !=============================================================================
33 : !=============================================================================
34 1536 : subroutine seasalt_init(seasalt_emis_scale)
35 : use sslt_sections, only: sslt_sections_init
36 : use constituents, only: cnst_get_ind
37 : use rad_constituents, only: rad_cnst_get_info
38 :
39 : real(r8), intent(in) :: seasalt_emis_scale
40 : integer :: m, l, nspec, ndx
41 : character(len=32) :: spec_name
42 :
43 1536 : seasalt_nbin = nslt
44 1536 : seasalt_nnum = nslt
45 4608 : allocate(seasalt_names(2*nslt))
46 4608 : allocate(seasalt_indices(2*nslt))
47 :
48 1536 : ndx=0
49 7680 : do m = 1, ntot_amode
50 6144 : call rad_cnst_get_info(0, m, nspec=nspec)
51 30720 : do l = 1, nspec
52 23040 : call rad_cnst_get_info(0, m, l, spec_name=spec_name )
53 29184 : if (spec_name(:3) == 'ncl') then
54 4608 : ndx=ndx+1
55 4608 : seasalt_names(ndx) = spec_name
56 4608 : seasalt_names(nslt+ndx) = 'num_'//spec_name(5:)
57 4608 : call cnst_get_ind(seasalt_names( ndx), seasalt_indices( ndx))
58 4608 : call cnst_get_ind(seasalt_names(nslt+ndx), seasalt_indices(nslt+ndx))
59 : endif
60 : enddo
61 : enddo
62 :
63 1536 : seasalt_active = any(seasalt_indices(:) > 0)
64 1536 : if (.not.seasalt_active) return
65 :
66 1536 : call sslt_sections_init()
67 :
68 1536 : emis_scale = seasalt_emis_scale
69 :
70 1536 : end subroutine seasalt_init
71 :
72 : !=============================================================================
73 : !=============================================================================
74 1489176 : subroutine seasalt_emis( u10cubed, srf_temp, ocnfrc, ncol, cflx )
75 :
76 1536 : use sslt_sections, only: nsections, fluxes, Dg, rdry
77 : use mo_constants, only: dns_aer_sst=>seasalt_density, pi
78 :
79 : ! dummy arguments
80 : real(r8), intent(in) :: u10cubed(:)
81 : real(r8), intent(in) :: srf_temp(:)
82 : real(r8), intent(in) :: ocnfrc(:)
83 : integer, intent(in) :: ncol
84 : real(r8), intent(inout) :: cflx(:,:)
85 :
86 : ! local vars
87 : integer :: mn, mm, ibin, isec, i
88 2978352 : real(r8) :: fi(ncol,nsections)
89 :
90 2978352 : real(r8) :: sst_sz_range_lo (nslt)
91 2978352 : real(r8) :: sst_sz_range_hi (nslt)
92 :
93 1489176 : if (nslt==4) then
94 0 : sst_sz_range_lo (:) = (/ 0.08e-6_r8, 0.02e-6_r8, 0.3e-6_r8, 1.0e-6_r8 /) ! accu, aitken, fine, coarse
95 0 : sst_sz_range_hi (:) = (/ 0.3e-6_r8, 0.08e-6_r8, 1.0e-6_r8, 10.0e-6_r8 /)
96 1489176 : else if (nslt==3) then
97 5956704 : sst_sz_range_lo (:) = (/ 0.08e-6_r8, 0.02e-6_r8, 1.0e-6_r8 /) ! accu, aitken, coarse
98 5956704 : sst_sz_range_hi (:) = (/ 1.0e-6_r8, 0.08e-6_r8, 10.0e-6_r8 /)
99 : endif
100 :
101 1489176 : fi(:ncol,:nsections) = fluxes( srf_temp, u10cubed, ncol )
102 :
103 5956704 : do ibin = 1,nslt
104 4467528 : mm = seasalt_indices(ibin)
105 4467528 : mn = seasalt_indices(nslt+ibin)
106 :
107 4467528 : if (mn>0) then
108 142960896 : do i=1, nsections
109 142960896 : if (Dg(i).ge.sst_sz_range_lo(ibin) .and. Dg(i).lt.sst_sz_range_hi(ibin)) then
110 671375952 : cflx(:ncol,mn)=cflx(:ncol,mn)+fi(:ncol,i)*ocnfrc(:ncol)*emis_scale !++ ag: scale sea-salt
111 : endif
112 : enddo
113 : endif
114 :
115 74597328 : cflx(:ncol,mm)=0.0_r8
116 144450072 : do i=1, nsections
117 142960896 : if (Dg(i).ge.sst_sz_range_lo(ibin) .and. Dg(i).lt.sst_sz_range_hi(ibin)) then
118 40207752 : cflx(:ncol,mm)=cflx(:ncol,mm)+fi(:ncol,i)*ocnfrc(:ncol)*emis_scale & !++ ag: scale sea-salt
119 711583704 : *4._r8/3._r8*pi*rdry(i)**3*dns_aer_sst ! should use dry size, convert from number to mass flux (kg/m2/s)
120 : endif
121 : enddo
122 :
123 : enddo
124 :
125 1489176 : end subroutine seasalt_emis
126 :
127 : end module seasalt_model
|