Line data Source code
1 : module hygrowghtpct_aerosol_optics_mod
2 :
3 : use shr_kind_mod, only: r8 => shr_kind_r8
4 : use aerosol_optics_mod, only: aerosol_optics
5 : use aerosol_state_mod, only: aerosol_state
6 : use aerosol_properties_mod, only: aerosol_properties
7 : use table_interp_mod, only: table_interp, table_interp_wghts, table_interp_calcwghts
8 :
9 : implicit none
10 :
11 : private
12 : public :: hygrowghtpct_aerosol_optics
13 :
14 : !> hygrowghtpct_aerosol_optics
15 : !! Table look up implementation of aerosol_optics to parameterize aerosol
16 : !! radiative properties in terms of weight precent of H2SO4/H2O solution
17 : type, extends(aerosol_optics) :: hygrowghtpct_aerosol_optics
18 :
19 : real(r8), allocatable :: totalmmr(:,:) ! total mmr of the aerosol
20 : real(r8), allocatable :: wgtpct(:,:) ! weight precent of H2SO4/H2O solution
21 :
22 : real(r8), pointer :: sw_hygro_ext_wtp(:,:) ! short wave extinction table
23 : real(r8), pointer :: sw_hygro_ssa_wtp(:,:) ! short wave single-scatter albedo table
24 : real(r8), pointer :: sw_hygro_asm_wtp(:,:) ! short wave asymmetry table
25 : real(r8), pointer :: lw_hygro_abs_wtp(:,:) ! long wave absorption table
26 :
27 : real(r8), pointer :: tbl_wgtpct(:) ! weight precent dimenstion values
28 :
29 : integer :: nwtp ! weight precent dimenstion size
30 :
31 : contains
32 :
33 : procedure :: sw_props
34 : procedure :: lw_props
35 :
36 : final :: destructor
37 :
38 : end type hygrowghtpct_aerosol_optics
39 :
40 : interface hygrowghtpct_aerosol_optics
41 : procedure :: constructor
42 : end interface hygrowghtpct_aerosol_optics
43 :
44 : contains
45 :
46 : !------------------------------------------------------------------------------
47 : !------------------------------------------------------------------------------
48 0 : function constructor(aero_props, aero_state, ilist, ibin, ncol, nlev, wgtpct_in) result(newobj)
49 :
50 : class(aerosol_properties),intent(in) :: aero_props ! aerosol_properties object
51 : class(aerosol_state),intent(in) :: aero_state ! aerosol_state object
52 : integer, intent(in) :: ilist ! climate or a diagnostic list number
53 : integer, intent(in) :: ibin ! bin number
54 : integer, intent(in) :: ncol ! number of columns
55 : integer, intent(in) :: nlev ! number of levels
56 : real(r8),intent(in) :: wgtpct_in(ncol,nlev) ! sulfate weight percent
57 :
58 : type(hygrowghtpct_aerosol_optics), pointer :: newobj
59 :
60 : integer :: ierr, nspec
61 : integer :: ispec
62 : integer :: i,k
63 :
64 0 : real(r8), pointer :: specmmr(:,:) ! species mass mixing ratio
65 :
66 0 : allocate(newobj, stat=ierr)
67 0 : if (ierr/=0) then
68 0 : nullify(newobj)
69 : return
70 : end if
71 :
72 0 : allocate(newobj%totalmmr(ncol,nlev),stat=ierr)
73 0 : if (ierr/=0) then
74 0 : nullify(newobj)
75 0 : return
76 : end if
77 :
78 0 : allocate(newobj%wgtpct(ncol,nlev),stat=ierr)
79 0 : if (ierr/=0) then
80 0 : nullify(newobj)
81 0 : return
82 : end if
83 :
84 : ! weight precent of H2SO4/H2O solution
85 0 : newobj%wgtpct(:ncol,:nlev) = wgtpct_in(:ncol,:nlev)
86 :
87 0 : call aero_props%optics_params(ilist, ibin, wgtpct=newobj%tbl_wgtpct, nwtp=newobj%nwtp)
88 :
89 0 : nspec = aero_props%nspecies(ilist, ibin)
90 :
91 0 : newobj%totalmmr(:,:) = 0._r8
92 :
93 0 : do ispec = 1,nspec
94 :
95 0 : call aero_state%get_ambient_mmr(ilist,ispec,ibin,specmmr)
96 0 : newobj%totalmmr(:ncol,:nlev) = newobj%totalmmr(:ncol,:nlev) + specmmr(:ncol,:nlev)
97 :
98 : end do
99 :
100 : call aero_props%optics_params(ilist, ibin, &
101 : sw_hygro_ext_wtp=newobj%sw_hygro_ext_wtp, &
102 : sw_hygro_ssa_wtp=newobj%sw_hygro_ssa_wtp, &
103 : sw_hygro_asm_wtp=newobj%sw_hygro_asm_wtp, &
104 0 : lw_hygro_ext_wtp=newobj%lw_hygro_abs_wtp)
105 :
106 0 : end function constructor
107 :
108 : !------------------------------------------------------------------------------
109 : ! returns short wave aerosol optics properties
110 : !------------------------------------------------------------------------------
111 0 : subroutine sw_props(self, ncol, ilev, iwav, pext, pabs, palb, pasm)
112 :
113 : class(hygrowghtpct_aerosol_optics), intent(in) :: self
114 : integer, intent(in) :: ncol ! number of columns
115 : integer, intent(in) :: ilev ! vertical level index
116 : integer, intent(in) :: iwav ! wave length index
117 : real(r8),intent(out) :: pext(ncol) ! parameterized specific extinction (m2/kg)
118 : real(r8),intent(out) :: pabs(ncol) ! parameterized specific absorption (m2/kg)
119 : real(r8),intent(out) :: palb(ncol) ! parameterized asymmetry factor
120 : real(r8),intent(out) :: pasm(ncol) ! parameterized single scattering albedo
121 :
122 : integer :: icol
123 0 : type(table_interp_wghts) :: wghts(ncol)
124 :
125 0 : wghts = table_interp_calcwghts( self%nwtp, self%tbl_wgtpct, ncol, self%wgtpct(:ncol,ilev) )
126 0 : pext = table_interp( ncol, self%nwtp, wghts, self%sw_hygro_ext_wtp(:,iwav) )
127 0 : pabs = (1._r8 - table_interp( ncol, self%nwtp, wghts, self%sw_hygro_ssa_wtp(:,iwav)))*pext
128 0 : pasm = table_interp( ncol, self%nwtp, wghts, self%sw_hygro_asm_wtp(:,iwav) )
129 :
130 0 : do icol = 1, ncol
131 :
132 0 : pext(icol) = pext(icol)*self%totalmmr(icol,ilev)
133 0 : pabs(icol) = pabs(icol)*self%totalmmr(icol,ilev)
134 0 : pabs(icol) = max(0._r8,pabs(icol))
135 0 : pabs(icol) = min(pext(icol),pabs(icol))
136 :
137 0 : palb(icol) = 1._r8-pabs(icol)/max(pext(icol),1.e-40_r8)
138 :
139 : end do
140 :
141 0 : end subroutine sw_props
142 :
143 : !------------------------------------------------------------------------------
144 : ! returns long wave aerosol optics properties
145 : !------------------------------------------------------------------------------
146 0 : subroutine lw_props(self, ncol, ilev, iwav, pabs)
147 :
148 : class(hygrowghtpct_aerosol_optics), intent(in) :: self
149 : integer, intent(in) :: ncol ! number of columns
150 : integer, intent(in) :: ilev ! vertical level index
151 : integer, intent(in) :: iwav ! wave length index
152 : real(r8),intent(out) :: pabs(ncol) ! parameterized specific absorption (m2/kg)
153 :
154 : integer :: icol
155 0 : type(table_interp_wghts) :: wghts(ncol)
156 :
157 0 : wghts = table_interp_calcwghts( self%nwtp, self%tbl_wgtpct, ncol, self%wgtpct(:ncol,ilev) )
158 :
159 0 : pabs = table_interp( ncol, self%nwtp, wghts, self%lw_hygro_abs_wtp(:,iwav) )
160 :
161 0 : do icol = 1, ncol
162 :
163 0 : pabs(icol) = pabs(icol)*self%totalmmr(icol,ilev)
164 0 : pabs(icol) = max(0._r8,pabs(icol))
165 :
166 : end do
167 :
168 0 : end subroutine lw_props
169 :
170 : !------------------------------------------------------------------------------
171 : !------------------------------------------------------------------------------
172 0 : subroutine destructor(self)
173 :
174 : type(hygrowghtpct_aerosol_optics), intent(inout) :: self
175 :
176 0 : deallocate(self%totalmmr)
177 0 : deallocate(self%wgtpct)
178 :
179 0 : nullify(self%tbl_wgtpct)
180 0 : nullify(self%sw_hygro_ext_wtp)
181 0 : nullify(self%sw_hygro_ssa_wtp)
182 0 : nullify(self%sw_hygro_asm_wtp)
183 0 : nullify(self%lw_hygro_abs_wtp)
184 :
185 0 : end subroutine destructor
186 :
187 0 : end module hygrowghtpct_aerosol_optics_mod
|