Line data Source code
1 : !-------------------------------------------------------------------------------
2 : ! This module uses the Lean solar irradiance data to provide a solar cycle
3 : ! scaling factor used in heating rate calculations
4 : !-------------------------------------------------------------------------------
5 : module rad_solar_var
6 :
7 : use shr_kind_mod , only : r8 => shr_kind_r8
8 : use solar_irrad_data, only : sol_irrad, we, nbins, do_spctrl_scaling, ssi_ref
9 : use solar_irrad_data, only : has_spectrum, has_ref_spectrum
10 : use cam_abortutils, only : endrun
11 :
12 : implicit none
13 : save
14 :
15 : private
16 : public :: rad_solar_var_init
17 : public :: get_variability
18 :
19 : real(r8), allocatable :: ref_band_irrad(:) ! scaling will be relative to ref_band_irrad in each band
20 : real(r8), allocatable :: irrad(:) ! solar irradiance at model timestep in each band
21 :
22 : real(r8), allocatable :: radbinmax(:)
23 : real(r8), allocatable :: radbinmin(:)
24 : integer :: nradbins
25 :
26 : !-------------------------------------------------------------------------------
27 : contains
28 : !-------------------------------------------------------------------------------
29 :
30 1536 : subroutine rad_solar_var_init( )
31 : use radconstants, only : get_number_sw_bands
32 : use radconstants, only : get_sw_spectral_boundaries
33 : use radconstants, only : get_ref_solar_band_irrad
34 : use radconstants, only : get_ref_total_solar_irrad
35 :
36 : integer :: i
37 : integer :: ierr
38 : integer :: yr, mon, tod
39 :
40 :
41 1536 : call get_number_sw_bands(nradbins)
42 :
43 1536 : if ( do_spctrl_scaling ) then
44 :
45 0 : if ( .not.has_spectrum ) then
46 0 : call endrun('rad_solar_var_init: solar input file must have irradiance spectrum')
47 : endif
48 :
49 0 : if ( .not.has_ref_spectrum ) then
50 0 : call endrun('rad_solar_var_init: solar input file must have reference irradiance spectrum')
51 : endif
52 :
53 0 : allocate (radbinmax(nradbins),stat=ierr)
54 0 : if (ierr /= 0) then
55 0 : call endrun('rad_solar_var_init: Error allocating space for radbinmax')
56 : end if
57 :
58 0 : allocate (radbinmin(nradbins),stat=ierr)
59 0 : if (ierr /= 0) then
60 0 : call endrun('rad_solar_var_init: Error allocating space for radbinmin')
61 : end if
62 :
63 0 : allocate (ref_band_irrad(nradbins), stat=ierr)
64 0 : if (ierr /= 0) then
65 0 : call endrun('rad_solar_var_init: Error allocating space for ref_band_irrad')
66 : end if
67 :
68 0 : allocate (irrad(nradbins), stat=ierr)
69 0 : if (ierr /= 0) then
70 0 : call endrun('rad_solar_var_init: Error allocating space for irrad')
71 : end if
72 :
73 0 : call get_sw_spectral_boundaries(radbinmin, radbinmax, 'nm')
74 :
75 0 : call integrate_spectrum( nbins, nradbins, we, radbinmin, radbinmax, ssi_ref, ref_band_irrad)
76 :
77 : endif
78 :
79 1536 : end subroutine rad_solar_var_init
80 :
81 : !-------------------------------------------------------------------------------
82 : !-------------------------------------------------------------------------------
83 0 : subroutine get_variability( sfac )
84 :
85 : real(r8), intent(out) :: sfac(nradbins) ! scaling factors for CAM heating
86 :
87 : integer :: yr, mon, day, tod
88 :
89 0 : if ( do_spctrl_scaling ) then
90 :
91 0 : call integrate_spectrum( nbins, nradbins, we, radbinmin, radbinmax, sol_irrad, irrad)
92 :
93 0 : sfac(:nradbins) = irrad(:nradbins)/ref_band_irrad(:nradbins)
94 :
95 : else
96 :
97 0 : sfac(:nradbins) = 1._r8
98 :
99 : endif
100 :
101 0 : end subroutine get_variability
102 :
103 : !-------------------------------------------------------------------------------
104 : ! private method.........
105 : !-------------------------------------------------------------------------------
106 :
107 0 : subroutine integrate_spectrum( nsrc, ntrg, src_x, min_trg, max_trg, src, trg )
108 :
109 : use mo_util, only : rebin
110 :
111 : implicit none
112 :
113 : !---------------------------------------------------------------
114 : ! ... dummy arguments
115 : !---------------------------------------------------------------
116 : integer, intent(in) :: nsrc ! dimension source array
117 : integer, intent(in) :: ntrg ! dimension target array
118 : real(r8), intent(in) :: src_x(nsrc+1) ! source coordinates
119 : real(r8), intent(in) :: max_trg(ntrg) ! target coordinates
120 : real(r8), intent(in) :: min_trg(ntrg) ! target coordinates
121 : real(r8), intent(in) :: src(nsrc) ! source array
122 : real(r8), intent(out) :: trg(ntrg) ! target array
123 :
124 : !---------------------------------------------------------------
125 : ! ... local variables
126 : !---------------------------------------------------------------
127 : real(r8) :: trg_x(2), targ(1) ! target coordinates
128 : integer :: i
129 :
130 0 : do i = 1, ntrg
131 :
132 0 : trg_x(1) = min_trg(i)
133 0 : trg_x(2) = max_trg(i)
134 :
135 0 : call rebin( nsrc, 1, src_x, trg_x, src(1:nsrc), targ(:) )
136 : ! W/m2/nm --> W/m2
137 0 : trg( i ) = targ(1)*(trg_x(2)-trg_x(1))
138 :
139 : enddo
140 :
141 :
142 0 : end subroutine integrate_spectrum
143 :
144 : end module rad_solar_var
|