Line data Source code
1 : module radconstants
2 :
3 : ! This module contains constants that are specific to the radiative transfer
4 : ! code used in the CAM4 model.
5 :
6 : use shr_kind_mod, only: r8 => shr_kind_r8
7 : use cam_abortutils, only: endrun
8 :
9 : implicit none
10 : private
11 :
12 : ! public routines
13 :
14 : public :: get_number_sw_bands
15 : public :: get_sw_spectral_boundaries
16 : public :: get_lw_spectral_boundaries
17 : public :: get_ref_solar_band_irrad
18 : public :: get_true_ref_solar_band_irrad
19 : public :: get_ref_total_solar_irrad
20 : public :: get_solar_band_fraction_irrad
21 : public :: radconstants_init
22 : public :: rad_gas_index
23 :
24 : ! SHORTWAVE DATA
25 :
26 : ! number of shorwave spectral intervals
27 : integer, parameter, public :: nswbands = 19
28 :
29 : integer, parameter, public :: idx_sw_diag = 8 ! index to sw visible band
30 :
31 : ! *** For interface consistency -- used only by modal_aero_optics and assumes use of RRTMG
32 : ! Need to provide a function interface to avoid this hack.
33 : integer, parameter, public :: idx_nir_diag = 999 ! index to sw near infrared (778-1240 nm) band
34 : integer, parameter, public :: idx_uv_diag = 999 ! index to sw uv (345-441 nm) band
35 : ! *** For interface consistency
36 :
37 :
38 : integer, parameter, public :: idx_lw_diag = 2 ! index to (H20 window) LW band
39 :
40 : ! LONGWAVE DATA
41 :
42 : ! number of lw bands
43 : integer, public, parameter :: nlwbands = 7
44 : ! Index of volc. abs., H2O non-window
45 : integer, public, parameter :: idx_LW_H2O_NONWND=1
46 : ! Index of volc. abs., H2O window
47 : integer, public, parameter :: idx_LW_H2O_WINDOW=2
48 : ! Index of volc. cnt. abs. 0500--0650 cm-1
49 : integer, public, parameter :: idx_LW_0500_0650=3
50 : ! Index of volc. cnt. abs. 0650--0800 cm-1
51 : integer, public, parameter :: idx_LW_0650_0800=4
52 : ! Index of volc. cnt. abs. 0800--1000 cm-1
53 : integer, public, parameter :: idx_LW_0800_1000=5
54 : ! Index of volc. cnt. abs. 1000--1200 cm-1
55 : integer, public, parameter :: idx_LW_1000_1200=6
56 : ! Index of volc. cnt. abs. 1200--2000 cm-1
57 : integer, public, parameter :: idx_LW_1200_2000=7
58 :
59 : ! GASES TREATED BY RADIATION (line spectrae)
60 :
61 : ! gasses required by radiation
62 : integer, public, parameter :: gasnamelength = 5
63 : integer, public, parameter :: nradgas = 8
64 : character(len=gasnamelength), public, parameter :: gaslist(nradgas) &
65 : = (/'H2O ','O3 ', 'O2 ', 'CO2 ', 'N2O ', 'CH4 ', 'CFC11', 'CFC12'/)
66 :
67 : ! what is the minimum mass mixing ratio that can be supported by radiation implementation?
68 : real(r8), public, parameter :: minmmr(nradgas) &
69 : = epsilon(1._r8)
70 :
71 : ! Solar and SW data for CAMRT
72 :
73 : ! Set index for cloud particle properties based on the wavelength,
74 : ! according to A. Slingo (1989) equations 1-3:
75 : ! Use index 1 (0.25 to 0.69 micrometers) for visible
76 : ! Use index 2 (0.69 - 1.19 micrometers) for near-infrared
77 : ! Use index 3 (1.19 to 2.38 micrometers) for near-infrared
78 : ! Use index 4 (2.38 to 4.00 micrometers) for near-infrared
79 : integer, public, parameter :: indxsl(nswbands) = &
80 : (/ 1, 1, 1, 1, 1, &
81 : 1, 1, 1, 1, 2, 3, &
82 : 3, 3, 3, 4, 4, &
83 : 4, 4, 4/)
84 :
85 : ! minimum wavelength of band in micrometers
86 : real(r8), parameter :: wavmin(nswbands) = &
87 : (/ 0.200_r8, 0.245_r8, 0.265_r8, 0.275_r8, 0.285_r8, &
88 : 0.295_r8, 0.305_r8, 0.350_r8, 0.640_r8, 0.700_r8, 0.700_r8, &
89 : 0.700_r8, 0.700_r8, 0.700_r8, 0.700_r8, 0.700_r8, &
90 : 2.630_r8, 4.160_r8, 4.160_r8/)
91 :
92 : ! maximum wavelength of band in micrometers
93 : real(r8), parameter :: wavmax(nswbands) = &
94 : (/ 0.245_r8, 0.265_r8, 0.275_r8, 0.285_r8, 0.295_r8, &
95 : 0.305_r8, 0.350_r8, 0.640_r8, 0.700_r8, 5.000_r8, 5.000_r8, &
96 : 5.000_r8, 5.000_r8, 5.000_r8, 5.000_r8, 5.000_r8, &
97 : 2.860_r8, 4.550_r8, 4.550_r8/)
98 :
99 : ! Fraction of solar flux in each stated spectral interval
100 : real(r8), public, parameter :: frcsol(nswbands) = &
101 : (/ .001488_r8, .001389_r8, .001290_r8, .001686_r8, .002877_r8, &
102 : .003869_r8, .026336_r8, .360739_r8, .065392_r8, .526861_r8, &
103 : .526861_r8, .526861_r8, .526861_r8, .526861_r8, .526861_r8, &
104 : .526861_r8, .006239_r8, .001834_r8, .001834_r8/)
105 :
106 : ! Weight of h2o in spectral interval
107 : real(r8), public, parameter :: ph2o(nswbands) = &
108 : (/ .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
109 : .000_r8, .000_r8, .000_r8, .000_r8, .505_r8, &
110 : .210_r8, .120_r8, .070_r8, .048_r8, .029_r8, &
111 : .018_r8, .000_r8, .000_r8, .000_r8/)
112 :
113 : ! Weight of co2 in spectral interval
114 : real(r8), public, parameter :: pco2(nswbands) = &
115 : (/ .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
116 : .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
117 : .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
118 : .000_r8, 1.000_r8, .640_r8, .360_r8/)
119 :
120 : ! Weight of o2 in spectral interval
121 : real(r8), public, parameter :: po2(nswbands) = &
122 : (/ .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
123 : .000_r8, .000_r8, .000_r8, 1.000_r8, 1.000_r8, &
124 : .000_r8, .000_r8, .000_r8, .000_r8, .000_r8, &
125 : .000_r8, .000_r8, .000_r8, .000_r8/)
126 :
127 : real(r8) :: solfrac_true(nswbands)
128 :
129 : ! Longwave spectral band limits (cm-1)
130 : real(r8), private, parameter :: wavenumber1_longwave(nlwbands) = &
131 : (/10._r8,350._r8,500._r8,650._r8,800._r8,1000._r8,1200._r8/)
132 :
133 : ! Longwave spectral band limits (cm-1)
134 : real(r8), private, parameter :: wavenumber2_longwave(nlwbands) = &
135 : (/350._r8,500._r8,650._r8,800._r8,1000._r8,1200._r8,2000._r8/)
136 :
137 : contains
138 :
139 :
140 : !------------------------------------------------------------------------------
141 1536 : subroutine get_number_sw_bands(number_of_bands)
142 : ! number of solar (shortwave) bands in the radiation code
143 : integer, intent(out) :: number_of_bands
144 :
145 1536 : number_of_bands = nswbands
146 :
147 1536 : end subroutine get_number_sw_bands
148 :
149 : !------------------------------------------------------------------------------
150 0 : subroutine get_lw_spectral_boundaries(low_boundaries, high_boundaries, units)
151 : ! provide spectral boundaries of each longwave band
152 :
153 : real(r8), intent(out) :: low_boundaries(nlwbands), high_boundaries(nlwbands)
154 : character(*), intent(in) :: units ! requested units
155 :
156 0 : select case (units)
157 : case ('inv_cm','cm^-1','cm-1')
158 0 : low_boundaries = wavenumber1_longwave
159 0 : high_boundaries = wavenumber2_longwave
160 : case('m','meter','meters')
161 0 : low_boundaries = 1.e-2_r8/wavenumber2_longwave
162 0 : high_boundaries = 1.e-2_r8/wavenumber1_longwave
163 : case('nm','nanometer','nanometers')
164 0 : low_boundaries = 1.e7_r8/wavenumber2_longwave
165 0 : high_boundaries = 1.e7_r8/wavenumber1_longwave
166 : case('um','micrometer','micrometers','micron','microns')
167 0 : low_boundaries = 1.e4_r8/wavenumber2_longwave
168 0 : high_boundaries = 1.e4_r8/wavenumber1_longwave
169 : case('cm','centimeter','centimeters')
170 0 : low_boundaries = 1._r8/wavenumber2_longwave
171 0 : high_boundaries = 1._r8/wavenumber1_longwave
172 : case default
173 0 : call endrun('get_lw_spectral_boundaries: spectral units not acceptable'//units)
174 : end select
175 :
176 0 : end subroutine get_lw_spectral_boundaries
177 :
178 : !------------------------------------------------------------------------------
179 7443877 : subroutine get_sw_spectral_boundaries(low_boundaries, high_boundaries, units)
180 : ! provide spectral boundaries of each shortwave band
181 :
182 : real(r8), intent(out) :: low_boundaries(nswbands), high_boundaries(nswbands)
183 : character(*), intent(in) :: units ! requested units
184 :
185 0 : select case (units)
186 : case ('inv_cm','cm^-1','cm-1')
187 0 : low_boundaries = 1.e4_r8/wavmax
188 0 : high_boundaries = 1.e4_r8/wavmin
189 : case('m','meter','meters')
190 0 : low_boundaries = 1.e-6_r8*wavmin
191 0 : high_boundaries = 1.e-6_r8*wavmax
192 : case('nm','nanometer','nanometers')
193 0 : low_boundaries = 1.e3_r8*wavmin
194 0 : high_boundaries = 1.e3_r8*wavmax
195 : case('um','micrometer','micrometers','micron','microns')
196 7443877 : low_boundaries = wavmin
197 7443877 : high_boundaries = wavmax
198 : case('cm','centimeter','centimeters')
199 0 : low_boundaries = 1.e-4_r8*wavmin
200 0 : high_boundaries = 1.e-4_r8*wavmax
201 : case default
202 7443877 : call endrun('get_sw_spectral_boundaries: spectral units not acceptable'//units)
203 : end select
204 :
205 7443877 : end subroutine get_sw_spectral_boundaries
206 :
207 : !------------------------------------------------------------------------------
208 0 : subroutine get_ref_solar_band_irrad( band_irrad )
209 :
210 : ! solar irradiance in each band (W/m^2)
211 : real(r8), intent(out) :: band_irrad(nswbands)
212 :
213 0 : band_irrad = frcsol
214 :
215 0 : end subroutine get_ref_solar_band_irrad
216 :
217 : !------------------------------------------------------------------------------
218 1536 : subroutine radconstants_init()
219 : ! The last bands are implemented as scalings to the solar flux
220 : ! so the corresponding actual flux applied to the heating
221 : ! is different from the solar in that band. These are the
222 : ! actual solar flux applied to each band
223 :
224 : integer :: ns
225 : real(r8):: psf(nswbands) ! scaled fractional solar spectrum in each band applied to unitary heating
226 :
227 30720 : do ns = 1, nswbands
228 29184 : psf(ns) = 1.0_r8
229 29184 : if(ph2o(ns)/=0._r8) psf(ns) = psf(ns)*ph2o(ns)
230 29184 : if(pco2(ns)/=0._r8) psf(ns) = psf(ns)*pco2(ns)
231 29184 : if(po2 (ns)/=0._r8) psf(ns) = psf(ns)*po2 (ns)
232 30720 : solfrac_true(ns) = frcsol(ns)*psf(ns)
233 : enddo
234 :
235 1536 : end subroutine radconstants_init
236 :
237 :
238 : !------------------------------------------------------------------------------
239 0 : subroutine get_true_ref_solar_band_irrad( solfrac_true_out )
240 :
241 : ! solar irradiance in each band (W/m^2)
242 :
243 : real(r8), intent(out) :: solfrac_true_out(nswbands)
244 :
245 0 : solfrac_true_out(:) = solfrac_true(:)
246 :
247 0 : end subroutine get_true_ref_solar_band_irrad
248 :
249 : !------------------------------------------------------------------------------
250 0 : subroutine get_ref_total_solar_irrad(tsi)
251 : ! provide Total Solar Irradiance assumed by radiation
252 :
253 : real(r8), intent(out) :: tsi
254 : real(r8) :: solfrac_true(nswbands)
255 :
256 0 : call get_true_ref_solar_band_irrad( solfrac_true )
257 0 : tsi = sum(solfrac_true)
258 :
259 0 : end subroutine get_ref_total_solar_irrad
260 :
261 : !------------------------------------------------------------------------------
262 0 : subroutine get_solar_band_fraction_irrad(fractional_irradiance)
263 : ! provide fractional solar irradiance in each band
264 :
265 : ! fraction of solar irradiance in each band
266 : real(r8), intent(out) :: fractional_irradiance(1:nswbands)
267 : real(r8) :: tsi ! total solar irradiance
268 :
269 0 : fractional_irradiance = frcsol
270 :
271 0 : end subroutine get_solar_band_fraction_irrad
272 :
273 : !------------------------------------------------------------------------------
274 7501512 : integer function rad_gas_index(gasname)
275 :
276 : ! return the index in the gaslist array of the specified gasname
277 :
278 : character(len=*),intent(in) :: gasname
279 : integer :: igas
280 :
281 7501512 : rad_gas_index = -1
282 30018384 : do igas = 1, nradgas
283 30018384 : if (trim(gaslist(igas)).eq.trim(gasname)) then
284 7501512 : rad_gas_index = igas
285 7501512 : return
286 : endif
287 : enddo
288 0 : call endrun ("rad_gas_index: can not find gas with name "//gasname)
289 0 : end function rad_gas_index
290 :
291 : end module radconstants
|