Merge pull request #22 from wirc-sjsu/develop-w21
[WRF-Fire-merge.git] / phys / module_cam_mp_radconstants.F
blob805b769b9d99cc681af55e6d90af28a9bec3e161
1 #define WRF_PORT
2 #define MODAL_AERO
3 ! Updated to CESM1.0.3 (CAM5.1.01) by Balwinder.Singh@pnnl.gov
5 module radconstants
7 ! This module contains constants that are specific to the radiative transfer
8 ! code used in the CAM3 model.
10 use shr_kind_mod,   only: r8 => shr_kind_r8
11 #ifndef WRF_PORT
12   use abortutils,     only: endrun
13 #else
14   use module_cam_support,    only: endrun
15 #endif
17 implicit none
18 private
20 ! public routines
22    public :: get_number_sw_bands
23    public :: get_sw_spectral_boundaries
24    public :: get_ref_solar_band_irrad
25    public :: get_true_ref_solar_band_irrad
26    public :: get_ref_total_solar_irrad
27    public :: get_solar_band_fraction_irrad
28    public :: radconstants_init
29    public :: rad_gas_index
31 ! optics files specify a type.  What length is it?
32 integer, parameter, public :: ot_length = 32
34 ! SHORTWAVE DATA
36 ! number of shorwave spectral intervals
37 integer, parameter, public :: nswbands = 19 
39 integer, parameter, public :: idx_sw_diag = 8 ! index to sw visible band
40 integer, parameter, public :: idx_lw_diag = 2 ! index to (H20 window) LW band
43 ! Number of evenly spaced intervals in rh
44 ! The globality of this mesh may not be necessary
45 ! Perhaps it could be specific to the aerosol
46 ! But it is difficult to see how refined it must be
47 ! for lookup.  This value was found to be sufficient
48 ! for Sulfate and probably necessary to resolve the
49 ! high variation near rh = 1.  Alternative methods
50 ! were found to be too slow.
51 ! Optimal approach would be for cam to specify size of aerosol
52 ! based on each aerosol's characteristics.  Radiation 
53 ! should know nothing about hygroscopic growth!
54 integer, parameter, public :: nrh = 1000  
56 ! LONGWAVE DATA
58 ! number of lw bands
59 integer, public, parameter  :: nlwbands = 7
60 ! Index of volc. abs., H2O non-window
61 integer, public, parameter :: idx_LW_H2O_NONWND=1
62 ! Index of volc. abs., H2O window
63 integer, public, parameter :: idx_LW_H2O_WINDOW=2
64 ! Index of volc. cnt. abs. 0500--0650 cm-1
65 integer, public, parameter :: idx_LW_0500_0650=3
66 ! Index of volc. cnt. abs. 0650--0800 cm-1
67 integer, public, parameter :: idx_LW_0650_0800=4
68 ! Index of volc. cnt. abs. 0800--1000 cm-1
69 integer, public, parameter :: idx_LW_0800_1000=5
70 ! Index of volc. cnt. abs. 1000--1200 cm-1
71 integer, public, parameter :: idx_LW_1000_1200=6
72 ! Index of volc. cnt. abs. 1200--2000 cm-1
73 integer, public, parameter :: idx_LW_1200_2000=7
75 ! GASES TREATED BY RADIATION (line spectrae)
77 ! gasses required by radiation
78 integer, public, parameter :: gasnamelength = 5
79 integer, public, parameter :: nradgas = 8
80 character(len=gasnamelength), public, parameter :: gaslist(nradgas) &
81    = (/'H2O  ','O3   ', 'O2   ', 'CO2  ', 'N2O  ', 'CH4  ', 'CFC11', 'CFC12'/)
83 ! what is the minimum mass mixing ratio that can be supported by radiation implementation?
84 real(r8), public, parameter :: minmmr(nradgas) &
85    = epsilon(1._r8)
87 ! Solar and SW data for CAMRT
89    ! minimum wavelength of band in micrometers
90    real(r8), public, parameter :: wavmin(nswbands) = &
91         (/   .200_r8,    .245_r8,    .265_r8,    .275_r8,    .285_r8, &
92              .295_r8,  .305_r8,    .350_r8,    .640_r8,    .700_r8,    .701_r8, &
93              .701_r8,  .701_r8,    .701_r8,    .702_r8,    .702_r8, &
94             2.630_r8, 4.160_r8,   4.160_r8/)
96    real(r8), public, parameter :: wavmin_true(nswbands) = &
97         (/   .200_r8,    .245_r8,    .265_r8,    .275_r8,    .285_r8, &
98              .295_r8,  .305_r8,    .350_r8,    .640_r8,    .700_r8,    .700_r8, &
99              .700_r8,  .700_r8,    .700_r8,    .700_r8,    .700_r8, &
100             2.630_r8, 4.160_r8,   4.160_r8/)
102    ! maximum wavelength of band in micrometers
103    real(r8), public, parameter :: wavmax(nswbands) = &
104         (/   .245_r8,  .265_r8,    .275_r8,    .285_r8,    .295_r8, &
105              .305_r8,  .350_r8,    .640_r8,    .700_r8,   5.000_r8,   5.000_r8, &
106             5.000_r8, 5.000_r8,   5.000_r8,   5.000_r8,   5.000_r8, &
107             2.860_r8, 4.550_r8,   4.550_r8/)
109    ! Fraction of solar flux in each stated spectral interval
110    real(r8), public, parameter :: frcsol(nswbands) = &
111      (/ .001488_r8, .001389_r8, .001290_r8, .001686_r8, .002877_r8, &
112         .003869_r8, .026336_r8, .360739_r8, .065392_r8, .526861_r8, &
113         .526861_r8, .526861_r8, .526861_r8, .526861_r8, .526861_r8, &
114         .526861_r8, .006239_r8, .001834_r8, .001834_r8/)
116    ! Weight of h2o in spectral interval
117    real(r8), public, parameter :: ph2o(nswbands) = &
118              (/    .000_r8,    .000_r8,    .000_r8,    .000_r8,    .000_r8, &
119         .000_r8,   .000_r8,    .000_r8,    .000_r8,    .505_r8,     &
120         .210_r8,   .120_r8,    .070_r8,    .048_r8,    .029_r8,     &
121         .018_r8,   .000_r8,    .000_r8,    .000_r8/)
123    ! Weight of co2 in spectral interval
124    real(r8), public, parameter :: pco2(nswbands) = &
125              (/    .000_r8,    .000_r8,    .000_r8,    .000_r8,    .000_r8, &
126         .000_r8,   .000_r8,    .000_r8,    .000_r8,    .000_r8,     &
127         .000_r8,   .000_r8,    .000_r8,    .000_r8,    .000_r8,     &
128         .000_r8,  1.000_r8,    .640_r8,    .360_r8/)
130    ! Weight of o2  in spectral interval
131    real(r8), public, parameter :: po2(nswbands) = &
132              (/    .000_r8,    .000_r8,    .000_r8,    .000_r8,    .000_r8, &
133         .000_r8,   .000_r8,    .000_r8,   1.000_r8,   1.000_r8,     &
134         .000_r8,   .000_r8,    .000_r8,    .000_r8,    .000_r8,     &
135         .000_r8,   .000_r8,    .000_r8,    .000_r8/)
137    real(r8) :: solfrac_true(nswbands)
139 contains
142 !------------------------------------------------------------------------------
143 subroutine get_number_sw_bands(number_of_bands)
144    ! number of solar (shortwave) bands in the radiation code
145    integer, intent(out) :: number_of_bands
147    number_of_bands = nswbands
149 end subroutine get_number_sw_bands
151 !------------------------------------------------------------------------------
152 subroutine get_sw_spectral_boundaries(low_boundaries, high_boundaries, units)
153    ! provide spectral boundaries of each shortwave band
155    real(r8), intent(out) :: low_boundaries(nswbands), high_boundaries(nswbands)
156    character(*), intent(in) :: units ! requested units
158    select case (units)
159    case ('inv_cm','cm^-1','cm-1')
160       low_boundaries = 1.e4_r8/wavmax
161       high_boundaries = 1.e4_r8/wavmin_true
162    case('m')
163       low_boundaries = 1.e-6_r8*wavmin_true
164       high_boundaries = 1.e-6_r8*wavmax
165    case('nm')
166       low_boundaries = 1.e3_r8*wavmin_true
167       high_boundaries = 1.e3_r8*wavmax
168    case('micrometer','micron','um')
169       low_boundaries = wavmin_true
170       high_boundaries = wavmax
171    case default
172       call endrun('rad_constants.F90: spectral units not acceptable'//units)
173    end select
175 end subroutine get_sw_spectral_boundaries
177 !------------------------------------------------------------------------------
178 subroutine get_ref_solar_band_irrad( band_irrad )
180    ! solar irradiance in each band (W/m^2)
181    real(r8), intent(out) :: band_irrad(nswbands)
183    band_irrad = frcsol
185 end subroutine get_ref_solar_band_irrad
187 !------------------------------------------------------------------------------
188 subroutine radconstants_init()
189 ! The last bands are implemented as scalings to the solar flux
190 ! so the corresponding actual flux applied to the heating
191 ! is different from the solar in that band.  These are the
192 ! actual solar flux applied to each band
194    integer :: ns
195    real(r8):: psf(nswbands)      !  scaled fractional solar spectrum in each band applied to unitary heating
197    do ns = 1, nswbands
198       psf(ns) = 1.0_r8
199       if(ph2o(ns)/=0._r8) psf(ns) = psf(ns)*ph2o(ns)
200       if(pco2(ns)/=0._r8) psf(ns) = psf(ns)*pco2(ns)
201       if(po2 (ns)/=0._r8) psf(ns) = psf(ns)*po2 (ns)
202       solfrac_true(ns)   = frcsol(ns)*psf(ns) 
203     enddo
205 end subroutine radconstants_init
208 !------------------------------------------------------------------------------
209 subroutine get_true_ref_solar_band_irrad( solfrac_true_out )
211    ! solar irradiance in each band (W/m^2)
213    real(r8), intent(out) :: solfrac_true_out(nswbands)
215    solfrac_true_out(:) = solfrac_true(:)
217 end subroutine get_true_ref_solar_band_irrad
219 !------------------------------------------------------------------------------
220 subroutine get_ref_total_solar_irrad(tsi)
221    ! provide Total Solar Irradiance assumed by radiation
223    real(r8), intent(out) :: tsi
224    real(r8) :: solfrac_true(nswbands)
226    call get_true_ref_solar_band_irrad( solfrac_true )
227    tsi = sum(solfrac_true)
229 end subroutine get_ref_total_solar_irrad
231 !------------------------------------------------------------------------------
232 subroutine get_solar_band_fraction_irrad(fractional_irradiance)
233    ! provide fractional solar irradiance in each band
235    ! fraction of solar irradiance in each band
236    real(r8), intent(out) :: fractional_irradiance(1:nswbands)
237    real(r8) :: tsi ! total solar irradiance
239    fractional_irradiance = frcsol
241 end subroutine get_solar_band_fraction_irrad
243 !------------------------------------------------------------------------------
244 integer function rad_gas_index(gasname)
246    ! return the index in the gaslist array of the specified gasname
248    character(len=*),intent(in) :: gasname
249    integer :: igas
251    rad_gas_index = -1
252    do igas = 1, nradgas
253       if (trim(gaslist(igas)).eq.trim(gasname)) then
254          rad_gas_index = igas
255          return
256       endif
257    enddo
258    call endrun ("rad_gas_index: can not find gas with name "//gasname)
259 end function rad_gas_index
261 end module radconstants