Line data Source code
1 : !-------------------------------------------------------------------------------
2 : ! solar variability parameters -- space weather indices
3 : !-------------------------------------------------------------------------------
4 : module solar_parms_data
5 :
6 : use shr_kind_mod, only : r8 => shr_kind_r8, shr_kind_cl
7 : use input_data_utils, only : time_coordinate
8 : use infnan, only : nan, assignment(=)
9 :
10 : implicit none
11 :
12 : private
13 : save
14 :
15 : ! public interface
16 :
17 : public :: solar_parms_init
18 : public :: solar_parms_advance
19 :
20 : logical, public :: solar_parms_on = .false.
21 :
22 : ! time-interpolated quantities
23 :
24 : real(r8), public, protected :: solar_parms_f107
25 : real(r8), public, protected :: solar_parms_f107a
26 : real(r8), public, protected :: solar_parms_f107p ! previous day
27 : real(r8), public, protected :: solar_parms_kp
28 : real(r8), public, protected :: solar_parms_ap
29 :
30 : ! private data
31 :
32 : real(r8), allocatable :: f107_in(:)
33 : real(r8), allocatable :: f107a_in(:)
34 : real(r8), allocatable :: kp_in(:)
35 : real(r8), allocatable :: ap_in(:)
36 :
37 : type(time_coordinate) :: time_coord_curr ! for current model time interpolation
38 : type(time_coordinate) :: time_coord_prev ! for previous day time interpolation
39 :
40 : contains
41 :
42 1536 : subroutine solar_parms_init(filepath, fixed, fixed_ymd, fixed_tod)
43 : !---------------------------------------------------------------
44 : ! ... initialize solar parmaters
45 : !---------------------------------------------------------------
46 :
47 : use ioFileMod
48 : use error_messages, only: alloc_err
49 : use cam_pio_utils, only: cam_pio_openfile
50 : use pio, only: file_desc_t, var_desc_t, pio_get_var, &
51 : pio_inq_varid, pio_closefile, pio_nowrite
52 :
53 : !---------------------------------------------------------------
54 : ! arguments
55 : !---------------------------------------------------------------
56 : character(len=*), intent(in) :: filepath
57 : logical, intent(in) :: fixed
58 : integer, intent(in) :: fixed_ymd
59 : integer, intent(in) :: fixed_tod
60 :
61 : !---------------------------------------------------------------
62 : ! ... local variables
63 : !---------------------------------------------------------------
64 : type(file_desc_t) :: ncid
65 : type(var_desc_t) :: varid
66 : integer :: astat
67 : character(len=shr_kind_cl) :: locfn
68 : integer :: ierr
69 :
70 1536 : solar_parms_f107 = nan
71 1536 : solar_parms_f107a = nan
72 1536 : solar_parms_f107p = nan
73 1536 : solar_parms_kp = nan
74 1536 : solar_parms_ap = nan
75 :
76 1536 : solar_parms_on = (filepath.ne.'NONE')
77 :
78 1536 : if (.not.solar_parms_on) return
79 :
80 : !-----------------------------------------------------------------------
81 : ! ... readin the solar parms dataset
82 : !-----------------------------------------------------------------------
83 :
84 0 : call getfil(filepath, locfn, 0)
85 0 : call cam_pio_openfile ( ncid, locfn, PIO_NOWRITE)
86 :
87 : call time_coord_prev%initialize( filepath, fixed=fixed, fixed_ymd=fixed_ymd, fixed_tod=fixed_tod, &
88 0 : force_time_interp=.true., try_dates=.true., delta_days=-1._r8 )
89 :
90 : call time_coord_curr%initialize( filepath, fixed=fixed, fixed_ymd=fixed_ymd, fixed_tod=fixed_tod, &
91 0 : force_time_interp=.true., try_dates=.true. )
92 :
93 : !---------------------------------------------------------------
94 : ! ... allocate and read solar parms
95 : !---------------------------------------------------------------
96 : allocate( f107_in(time_coord_curr%ntimes), f107a_in(time_coord_curr%ntimes), &
97 0 : kp_in(time_coord_curr%ntimes), ap_in(time_coord_curr%ntimes), stat=astat )
98 0 : if( astat /= 0 ) then
99 0 : call alloc_err( astat, 'solar_parms_init', 'f107_in ... ap_in ', time_coord_curr%ntimes )
100 : end if
101 0 : ierr = pio_inq_varid( ncid, 'f107', varid )
102 0 : ierr = pio_get_var( ncid, varid, f107_in )
103 0 : ierr = pio_inq_varid( ncid, 'f107a', varid )
104 0 : ierr = pio_get_var( ncid, varid, f107a_in )
105 0 : ierr = pio_inq_varid( ncid, 'kp', varid )
106 0 : ierr = pio_get_var( ncid, varid, kp_in )
107 0 : ierr = pio_inq_varid( ncid, 'ap', varid )
108 0 : ierr = pio_get_var( ncid, varid, ap_in )
109 :
110 0 : call pio_closefile( ncid )
111 :
112 3072 : end subroutine solar_parms_init
113 :
114 370944 : subroutine solar_parms_advance
115 : !---------------------------------------------------------------
116 : ! time interpolate space wx indices
117 : !---------------------------------------------------------------
118 :
119 : integer :: ndx1, ndx2
120 : real(r8) :: wgt1, wgt2
121 :
122 370944 : if (solar_parms_on) then
123 0 : call time_coord_curr%advance()
124 0 : ndx1=time_coord_curr%indxs(1)
125 0 : ndx2=time_coord_curr%indxs(2)
126 0 : wgt1=time_coord_curr%wghts(1)
127 0 : wgt2=time_coord_curr%wghts(2)
128 0 : solar_parms_f107 = wgt1*f107_in(ndx1) + wgt2*f107_in(ndx2)
129 0 : solar_parms_f107a = wgt1*f107a_in(ndx1) + wgt2*f107a_in(ndx2)
130 0 : solar_parms_kp = wgt1*kp_in(ndx1) + wgt2*kp_in(ndx2)
131 0 : solar_parms_ap = wgt1*ap_in(ndx1) + wgt2*ap_in(ndx2)
132 :
133 0 : call time_coord_prev%advance()
134 0 : ndx1=time_coord_prev%indxs(1)
135 0 : ndx2=time_coord_prev%indxs(2)
136 0 : wgt1=time_coord_prev%wghts(1)
137 0 : wgt2=time_coord_prev%wghts(2)
138 0 : solar_parms_f107p = wgt1*f107_in(ndx1) + wgt2*f107_in(ndx2)
139 : endif
140 :
141 1536 : end subroutine solar_parms_advance
142 :
143 : end module solar_parms_data
|