Line data Source code
1 : !-------------------------------------------------------------------
2 : ! rebins the 4 sea salt bins into 2 bins for the radiation
3 : !
4 : ! N.B. This code looks for the constituents of SSLTA and SSLTC
5 : ! in the physics buffer first, and uses those if found.
6 : ! Consequently, it is not possible to have prognostic sea
7 : ! salt be radiatively active if the prescribed sea salt is
8 : ! also present. The current (cam3_5_52) chemistry configurations
9 : ! don't allow both prescribed and prognostic to be present
10 : ! simultaneously, but a more flexible chemistry package that
11 : ! allows this would break this code.
12 : !
13 : ! Created by: Francis Vitt
14 : ! Date: 9 May 2008
15 : !-------------------------------------------------------------------
16 : module sslt_rebin
17 :
18 : use shr_kind_mod, only: r8 => shr_kind_r8
19 :
20 : implicit none
21 :
22 : integer :: indices(4)
23 : integer :: sslta_idx, ssltc_idx
24 :
25 : logical :: has_sslt = .false.
26 : character(len=1) :: source
27 : character(len=1), parameter :: DATA = 'D'
28 : character(len=1), parameter :: PROG = 'P'
29 :
30 : private
31 : public :: sslt_rebin_init, sslt_rebin_adv, sslt_rebin_register
32 : contains
33 :
34 :
35 : !-------------------------------------------------------------------
36 : !-------------------------------------------------------------------
37 1496904 : subroutine sslt_rebin_register
38 : use ppgrid, only : pver,pcols
39 :
40 : use physics_buffer, only : pbuf_add_field, dtype_r8
41 :
42 : ! add SSLTA and SSLTC to physics buffer
43 1536 : call pbuf_add_field('SSLTA','physpkg',dtype_r8,(/pcols,pver/),sslta_idx)
44 1536 : call pbuf_add_field('SSLTC','physpkg',dtype_r8,(/pcols,pver/),ssltc_idx)
45 :
46 1536 : endsubroutine sslt_rebin_register
47 :
48 : !-------------------------------------------------------------------
49 : !-------------------------------------------------------------------
50 1536 : subroutine sslt_rebin_init()
51 :
52 1536 : use constituents, only : cnst_get_ind
53 :
54 : use physics_buffer, only : pbuf_get_index
55 : use ppgrid, only : pver
56 : use cam_history, only : addfld
57 :
58 : implicit none
59 :
60 : integer :: errcode
61 :
62 :
63 1536 : indices(1) = pbuf_get_index('sslt1',errcode)
64 1536 : indices(2) = pbuf_get_index('sslt2',errcode)
65 1536 : indices(3) = pbuf_get_index('sslt3',errcode)
66 1536 : indices(4) = pbuf_get_index('sslt4',errcode)
67 :
68 1536 : has_sslt = all( indices(:) > 0 )
69 1536 : if ( has_sslt ) source = DATA
70 :
71 1536 : if ( .not. has_sslt ) then
72 1536 : call cnst_get_ind ('SSLT01', indices(1), abort=.false.)
73 1536 : call cnst_get_ind ('SSLT02', indices(2), abort=.false.)
74 1536 : call cnst_get_ind ('SSLT03', indices(3), abort=.false.)
75 1536 : call cnst_get_ind ('SSLT04', indices(4), abort=.false.)
76 1536 : has_sslt = all( indices(:) > 0 )
77 1536 : if ( has_sslt ) source = PROG
78 : endif
79 :
80 1536 : if ( has_sslt ) then
81 0 : call addfld('SSLTA', (/ 'lev' /), 'A','kg/kg', 'sea salt' )
82 0 : call addfld('SSLTC', (/ 'lev' /), 'A','kg/kg', 'sea salt' )
83 : endif
84 :
85 1536 : end subroutine sslt_rebin_init
86 :
87 : !-------------------------------------------------------------------
88 : !-------------------------------------------------------------------
89 1495368 : subroutine sslt_rebin_adv(pbuf, phys_state)
90 :
91 1536 : use physics_types,only : physics_state
92 :
93 : use ppgrid, only : pver, pcols
94 : use cam_history, only : outfld
95 : use physics_buffer, only : physics_buffer_desc, pbuf_get_field
96 :
97 : implicit none
98 :
99 :
100 : type(physics_state), target, intent(in) :: phys_state
101 : type(physics_buffer_desc), pointer :: pbuf(:)
102 :
103 : !++ changed wgt_sscm declaration for roundoff validation with earlier code
104 : ! real(r8), parameter :: wgt_sscm = 6.0_r8 / 7.0_r8 ! Fraction of total seasalt mass in coarse mode
105 : real(r8), parameter :: wgt_sscm = 6.0_r8 / 7.0_r8 ! Fraction of total seasalt mass in coarse mode
106 :
107 1495368 : real(r8), dimension(:,:), pointer :: sslt1, sslt2, sslt3, sslt4
108 1495368 : real(r8), dimension(:,:), pointer :: sslta, ssltc
109 : integer :: lchnk, ncol
110 : real(r8) :: sslt_sum(pcols,pver)
111 :
112 1495368 : lchnk = phys_state%lchnk
113 1495368 : ncol = phys_state%ncol
114 :
115 1495368 : if (.not. has_sslt) return
116 :
117 0 : select case( source )
118 : case (PROG)
119 0 : sslt1 => phys_state%q(:,:,indices(1))
120 0 : sslt2 => phys_state%q(:,:,indices(2))
121 0 : sslt3 => phys_state%q(:,:,indices(3))
122 0 : sslt4 => phys_state%q(:,:,indices(4))
123 : case (DATA)
124 0 : call pbuf_get_field(pbuf, indices(1), sslt1)
125 0 : call pbuf_get_field(pbuf, indices(2), sslt2)
126 0 : call pbuf_get_field(pbuf, indices(3), sslt3)
127 0 : call pbuf_get_field(pbuf, indices(4), sslt4)
128 : end select
129 :
130 0 : call pbuf_get_field(pbuf, sslta_idx, sslta )
131 0 : call pbuf_get_field(pbuf, ssltc_idx, ssltc )
132 :
133 0 : sslt_sum(:ncol,:) = sslt1(:ncol,:) + sslt2(:ncol,:) + sslt3(:ncol,:) + sslt4(:ncol,:)
134 0 : sslta(:ncol,:) = (1._r8-wgt_sscm)*sslt_sum(:ncol,:) ! fraction of seasalt mass in accumulation mode
135 0 : ssltc(:ncol,:) = wgt_sscm*sslt_sum(:ncol,:) ! fraction of seasalt mass in coagulation mode
136 :
137 0 : call outfld( 'SSLTA', sslta(:ncol,:), ncol, lchnk )
138 0 : call outfld( 'SSLTC', ssltc(:ncol,:), ncol, lchnk )
139 :
140 1495368 : end subroutine sslt_rebin_adv
141 :
142 : end module sslt_rebin
|