IODA Bundle
define_mod.f90
Go to the documentation of this file.
1 module define_mod
2 
3 use kinds, only: r_kind, i_kind
5 use netcdf, only: nf90_float, nf90_int, nf90_char
6 
7 implicit none
8 
9 real(r_kind), parameter :: t_kelvin = 273.15
10 real(r_kind), parameter :: missing_r = -999.0
11 integer(i_kind), parameter :: missing_i = -999
12 integer(i_kind), parameter :: not_use = 100
13 integer(i_kind), parameter :: itrue = 1
14 integer(i_kind), parameter :: ifalse = 0
15 integer(i_kind), parameter :: nstring = 50
16 integer(i_kind), parameter :: ndatetime = 20
17 integer(i_kind), parameter :: nobtype = 6 ! number of ob types
18 integer(i_kind), parameter :: n_ncdim = 5 ! total numner of nc dimensions
19 integer(i_kind), parameter :: nvar_met = 6
20 integer(i_kind), parameter :: nvar_info = 9 ! number of metadata
21 integer(i_kind), parameter :: nsen_info = 7 ! number of sensor metadata
22 integer(i_kind), parameter :: ninst = 12
23 !integer(i_kind), parameter :: ninst = 13 ! including airs
24 integer(i_kind), parameter :: write_nc_conv = 1
25 integer(i_kind), parameter :: write_nc_radiance = 2
26 
27 ! variables for defining observation types and met variables each type has
28 character(len=nstring), dimension(nobtype) :: obtype_list = &
29  (/ &
30  'sondes ', &
31  'aircraft ', &
32  'sfc ', &
33  'satwind ', &
34  'profiler ', &
35  'ascat ' &
36  /)
37 
38 character(len=nstring), dimension(nvar_met) :: name_var_met = &
39  (/ &
40  var_u, &
41  var_v, &
42  var_ts, &
43  var_tv, &
44  var_q, &
45  var_ps &
46  /)
47 
48 ! variable flags for var_u, var_v, var_ts, var_tv, var_q, var_ps
49 integer(i_kind), dimension(nvar_met,nobtype) :: vflag = reshape ( &
50  (/ &
51  itrue, itrue, itrue, itrue, itrue, ifalse, & ! sonde
52  itrue, itrue, itrue, itrue, itrue, ifalse, & ! aircraft
53  itrue, itrue, itrue, itrue, itrue, itrue, & ! sfc
54  itrue, itrue, ifalse, ifalse, ifalse, ifalse, & ! satwnd
55  itrue, itrue, ifalse, ifalse, ifalse, ifalse, & ! profiler
56  itrue, itrue, ifalse, ifalse, ifalse, ifalse & ! ascat
57  /), (/nvar_met,nobtype/) )
58 
59 character(len=nstring), dimension(nvar_met) :: unit_var_met = &
60  (/ &
61  'm/s ', &
62  'm/s ', &
63  'K ', &
64  'K ', &
65  'kg/kg ', &
66  'Pa ' &
67  /)
68 
69 ! variables for defining radiance instrument types
70 character(len=nstring), dimension(ninst) :: inst_list = &
71  (/ &
72  'amsua_n15 ', &
73  'amsua_n18 ', &
74  'amsua_n19 ', &
75  'amsua_metop-a ', &
76  'amsua_metop-b ', &
77  'amsua_metop-c ', &
78 ! 'airs_aqua ', &
79  'amsua_aqua ', &
80  'mhs_n18 ', &
81  'mhs_n19 ', &
82  'mhs_metop-a ', &
83  'mhs_metop-b ', &
84  'mhs_metop-c ' &
85  /)
86 
87 ! variables for outputing netcdf files
88 character(len=nstring), dimension(n_ncdim) :: name_ncdim = &
89  (/ &
90  'nvars ', &
91  'nlocs ', &
92  'nrec ', &
93  'nstring ', &
94  'ndatetime ' &
95  /)
96 character(len=nstring), dimension(nvar_info) :: name_var_info = &
97  (/ &
98  'air_pressure ', &
99  'height ', &
100  'station_elevation', &
101  'latitude ', &
102  'longitude ', &
103  'record_number ', &
104  'datetime ', &
105  'station_id ', &
106  'variable_names ' &
107  /)
108 integer(i_kind), dimension(nvar_info) :: type_var_info = &
109  (/ &
110  nf90_float, &
111  nf90_float, &
112  nf90_float, &
113  nf90_float, &
114  nf90_float, &
115  nf90_int, &
116  nf90_char, &
117  nf90_char, &
118  nf90_char &
119  /)
120 character(len=nstring), dimension(2,nvar_info) :: dim_var_info = reshape ( &
121  (/ &
122  'nlocs ', 'null ', &
123  'nlocs ', 'null ', &
124  'nlocs ', 'null ', &
125  'nlocs ', 'null ', &
126  'nlocs ', 'null ', &
127  'nlocs ', 'null ', &
128  'ndatetime ', 'nlocs ', &
129  'nstring ', 'nlocs ', &
130  'nstring ', 'nvars ' &
131  /), (/2, nvar_info/) )
132 character(len=nstring), dimension(nsen_info) :: name_sen_info = &
133  (/ &
134  'solar_azimuth_angle ', &
135  'scan_position ', &
136  'sensor_azimuth_angle', &
137  'solar_zenith_angle ', &
138  'sensor_zenith_angle ', &
139  'sensor_view_angle ', &
140  'sensor_channel ' &
141  /)
142 integer(i_kind), dimension(nsen_info) :: type_sen_info = &
143  (/ &
144  nf90_float, &
145  nf90_float, &
146  nf90_float, &
147  nf90_float, &
148  nf90_float, &
149  nf90_float, &
150  nf90_int &
151  /)
152 character(len=nstring), dimension(2,nsen_info) :: dim_sen_info = reshape ( &
153  (/ &
154  'nlocs ', 'null ', &
155  'nlocs ', 'null ', &
156  'nlocs ', 'null ', &
157  'nlocs ', 'null ', &
158  'nlocs ', 'null ', &
159  'nlocs ', 'null ', &
160  'nvars ', 'null ' &
161  /), (/2,nsen_info/) )
162 
163 ! variables for storing data
165  real(r_kind) :: val ! observation value
166  integer(i_kind) :: qm ! observation quality marker
167  real(r_kind) :: err ! observational error
168  integer(i_kind) :: rptype ! report type
169 end type xfield_type
170 
172  integer(i_kind) :: nvars
173  integer(i_kind) :: nrecs
174  integer(i_kind) :: nlocs
175  integer(i_kind), allocatable, dimension(:) :: var_idx
176  type (xfield_type), allocatable, dimension(:,:) :: xfield
177  real(r_kind), allocatable, dimension(:,:) :: xinfo_float
178  integer(i_kind), allocatable, dimension(:,:) :: xinfo_int
179  character(len=nstring), allocatable, dimension(:,:) :: xinfo_char
180  real(r_kind), allocatable, dimension(:,:) :: xseninfo_float
181  integer(i_kind), allocatable, dimension(:,:) :: xseninfo_int
182  character(len=nstring), allocatable, dimension(:,:) :: xseninfo_char
183 end type xdata_type
184 
185 type(xdata_type), allocatable, dimension(:) :: xdata ! dim: number of ob types
186 
187 contains
188 
189 subroutine set_obtype_conv(t29, obtype)
190 
191 ! https://www.emc.ncep.noaa.gov/BUFRLIB/tables/CodeFlag_0_STDv33_LOC7.html#055008
192 
193 ! assign conventional obtype name based on data dump report type (t29)
194 ! obtype names here should be consistent with those defined in obtype_list
195 
196  implicit none
197 
198  integer(i_kind), intent(in) :: t29
199  character(len=*), intent(out) :: obtype
200 
201  obtype = 'unknown'
202 
203  select case(t29)
204  case (11, 12, 13, 22, 23, 31)
205  obtype = 'sondes'
206  !select case (kx)
207  !case (120, 122, 132, 220, 222, 232)
208  ! obtype = 'sondes'
209  !case (221)
210  ! obtype = 'pilot'
211  !end select
212  case (41)
213  ! kx case (130:131, 133, 230:231, 233)
214  obtype = 'aircraft'
215  case (522, 523)
216  obtype = 'sfc' ! 'ship'
217  case (531, 532, 561, 562, 534)
218  obtype = 'sfc' ! 'buoy'
219  case (511, 514, 540) ! mesonet 540
220  ! kx case (181, 281)
221  obtype = 'sfc' ! 'synop'
222  case (512)
223  ! kx case (187, 287)
224  obtype = 'sfc' ! 'metar'
225  case (63)
226  ! kx case (242:246, 252:253, 255)
227  obtype = 'satwind'
228  case (581, 582, 583, 584)
229  ! ERS 581, QuikSCAT 582, WindSat 583, ASCAT 584
230  obtype = 'ascat'
231  !case (74)
232  ! obtype = 'gpspw'
233  !case (71, 73, 75, 76, 77)
234  ! t29=73/kx=229 (Wind profiler originating in PIBAL bulletins (tropical and European)
235  ! t29=77/kx=126 (Multi-Agency Profiler (MAP) RASS temperatures)
236  case (71, 73, 75, 76)
237  obtype = 'profiler'
238  !case (571, 65)
239  ! obtype = 'ssmir' ! ssmi retrieval
240  end select
241 
242 end subroutine set_obtype_conv
243 
244 subroutine set_name_satellite(satid, satellite)
245 
246 ! https://www.emc.ncep.noaa.gov/BUFRLIB/tables/CodeFlag_0_STDv33_LOC7.html#001007
247 
248 ! assign satellite name based on BUFR SAID (0-01-007)
249 
250  implicit none
251 
252  integer(i_kind), intent(in) :: satid
253  character(len=*), intent(out) :: satellite
254 
255  satellite = 'unknown'
256 
257  select case ( satid )
258  case ( 3 ); satellite = 'metop-b'
259  case ( 4 ); satellite = 'metop-a'
260  case ( 5 ); satellite = 'metop-c'
261  case ( 206 ); satellite = 'n15'
262  case ( 207 ); satellite = 'n16'
263  case ( 208 ); satellite = 'n17'
264  case ( 209 ); satellite = 'n18'
265  case ( 223 ); satellite = 'n19'
266  case ( 224 ); satellite = 'npp'
267  case ( 225 ); satellite = 'n20'
268  case ( 226 ); satellite = 'n21'
269  case ( 783 ); satellite = 'terra'
270  case ( 784 ); satellite = 'aqua'
271  end select
272 
273 end subroutine set_name_satellite
274 
275 subroutine set_name_sensor(instid, sensor)
276 
277 ! https://www.emc.ncep.noaa.gov/BUFRLIB/tables/CodeFlag_0_STDv33_LOC7.html#002019
278 
279 ! assign sensor name based on BUFR SIID (0-02-019) Satellite instruments
280 
281  implicit none
282 
283  integer(i_kind), intent(in) :: instid
284  character(len=*), intent(out) :: sensor
285 
286  sensor = 'unknown'
287 
288  select case ( instid )
289  case ( 570 ); sensor = 'amsua'
290  case ( 574 ); sensor = 'amsub'
291  case ( 203 ); sensor = 'mhs'
292  case ( 617 ); sensor = 'abi'
293  case ( 297 ); sensor = 'ahi'
294  case ( 207 ); sensor = 'seviri'
295  case ( 420 ); sensor = 'airs'
296  case ( 620 ); sensor = 'cris'
297  case ( 221 ); sensor = 'iasi'
298  case ( 389 ); sensor = 'modis'
299  case ( 616 ); sensor = 'viirs'
300  end select
301 
302 end subroutine set_name_sensor
303 
304 subroutine set_brit_obserr(name_inst, ichan, obserr)
305 
306 ! set brightness temperature observation errors
307 
308 ! For now it is a temporary subroutine to assign AMSU-A and MHS observation errors
309 
310  implicit none
311 
312  character(len=*), intent(in) :: name_inst ! instrument name eg. amsua_n15
313  integer(i_kind), intent(in) :: ichan ! channel index
314  real(r_kind), intent(out) :: obserr ! obserr of ichan
315 
316  integer(i_kind) :: nchan
317  real(r_kind), allocatable, dimension(:) :: obserrors
318 
319  obserr = missing_r
320 
321  if ( name_inst(1:5) == 'amsua' ) then
322  nchan = 15
323  allocate(obserrors(nchan))
324  select case ( trim(name_inst) )
325  case ( 'amsua_n15' )
326  obserrors = (/ 3.0, 2.2, 2.0, 0.6, 0.3, 0.23, 0.25, 0.275, 0.34, 0.4, 0.6, 1.0, 1.5, 2.0, 3.5 /)
327  case ( 'amsua_n18' )
328  obserrors = (/ 2.5, 2.2, 2.0, 0.55, 0.3, 0.23, 0.23, 0.25, 0.25, 0.35, 0.4, 0.55, 0.8, 3.0, 3.5 /)
329  case ( 'amsua_n19' )
330  obserrors = (/ 2.5, 2.2, 2.0, 0.55, 0.3, 0.23, 0.23, 0.25, 0.25, 0.35, 0.4, 0.55, 0.8, 3.0, 3.5 /)
331  case ( 'amsua_metop-b' )
332  obserrors = (/ 2.5, 2.2, 2.0, 0.55, 0.3, 0.23, 0.23, 0.25, 0.25, 0.35, 0.4, 0.55, 0.8, 3.0, 3.5 /)
333  case ( 'amsua_metop-a' )
334  obserrors = (/ 2.5, 2.2, 2.0, 0.55, 0.3, 0.23, 0.23, 0.25, 0.25, 0.35, 0.4, 0.55, 0.8, 3.0, 3.5 /)
335  case ( 'amsua_metop-c' )
336  obserrors = (/ 2.5, 2.2, 2.0, 0.55, 0.3, 0.23, 0.23, 0.25, 0.25, 0.35, 0.4, 0.55, 0.8, 3.0, 3.5 /)
337  case ( 'amsua_aqua' )
338  obserrors = (/ 2.5, 2.0, 2.0, 0.5, 0.4, 0.4, 0.5, 0.3, 0.35, 0.35, 0.45, 1.0, 1.5, 2.5, 2.5 /)
339  case default
340  return
341  end select
342  else if ( name_inst(1:3) == 'mhs' ) then
343  nchan = 5
344  allocate(obserrors(nchan))
345  select case ( trim(name_inst) )
346  case ( 'mhs_n18' )
347  obserrors = (/ 2.5, 2.5, 2.5, 2.0, 2.0 /)
348  case ( 'mhs_n19' )
349  obserrors = (/ 2.5, 2.5, 2.5, 2.0, 2.0 /)
350  case ( 'mhs_metop-a' )
351  obserrors = (/ 2.5, 2.5, 2.5, 2.0, 2.0 /)
352  case ( 'mhs_metop-b' )
353  obserrors = (/ 2.5, 2.5, 2.5, 2.0, 2.0 /)
354  case ( 'mhs_metop-c' )
355  obserrors = (/ 2.5, 2.5, 2.5, 2.0, 2.0 /)
356  case default
357  return
358  end select
359  else
360  return
361  end if
362 
363  obserr = obserrors(ichan)
364 
365  deallocate(obserrors)
366 
367 end subroutine set_brit_obserr
368 
369 end module define_mod
integer(i_kind), parameter nstring
Definition: define_mod.f90:15
integer(i_kind), parameter nvar_met
Definition: define_mod.f90:19
integer(i_kind), parameter write_nc_radiance
Definition: define_mod.f90:25
character(len=nstring), dimension(nobtype) obtype_list
Definition: define_mod.f90:28
character(len=nstring), dimension(2, nsen_info) dim_sen_info
Definition: define_mod.f90:152
character(len=nstring), dimension(nvar_info) name_var_info
Definition: define_mod.f90:96
integer(i_kind), parameter nobtype
Definition: define_mod.f90:17
integer(i_kind), dimension(nsen_info) type_sen_info
Definition: define_mod.f90:142
integer(i_kind), parameter itrue
Definition: define_mod.f90:13
character(len=nstring), dimension(n_ncdim) name_ncdim
Definition: define_mod.f90:88
subroutine set_name_sensor(instid, sensor)
Definition: define_mod.f90:276
real(r_kind), parameter t_kelvin
Definition: define_mod.f90:9
integer(i_kind), parameter ifalse
Definition: define_mod.f90:14
integer(i_kind), parameter write_nc_conv
Definition: define_mod.f90:24
integer(i_kind), parameter not_use
Definition: define_mod.f90:12
integer(i_kind), parameter ninst
Definition: define_mod.f90:22
integer(i_kind), parameter nsen_info
Definition: define_mod.f90:21
integer(i_kind), parameter missing_i
Definition: define_mod.f90:11
integer(i_kind), parameter ndatetime
Definition: define_mod.f90:16
character(len=nstring), dimension(nsen_info) name_sen_info
Definition: define_mod.f90:132
integer(i_kind), parameter nvar_info
Definition: define_mod.f90:20
character(len=nstring), dimension(nvar_met) name_var_met
Definition: define_mod.f90:38
subroutine set_obtype_conv(t29, obtype)
Definition: define_mod.f90:190
real(r_kind), parameter missing_r
Definition: define_mod.f90:10
character(len=nstring), dimension(ninst) inst_list
Definition: define_mod.f90:70
type(xdata_type), dimension(:), allocatable xdata
Definition: define_mod.f90:185
integer(i_kind), dimension(nvar_met, nobtype) vflag
Definition: define_mod.f90:49
subroutine set_name_satellite(satid, satellite)
Definition: define_mod.f90:245
integer(i_kind), parameter n_ncdim
Definition: define_mod.f90:18
character(len=nstring), dimension(nvar_met) unit_var_met
Definition: define_mod.f90:59
subroutine set_brit_obserr(name_inst, ichan, obserr)
Definition: define_mod.f90:305
integer(i_kind), dimension(nvar_info) type_var_info
Definition: define_mod.f90:108
character(len=nstring), dimension(2, nvar_info) dim_var_info
Definition: define_mod.f90:120
Definition: kinds.f90:1
integer, parameter, public i_kind
Definition: kinds.f90:71
integer, parameter, public r_kind
Definition: kinds.f90:100
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_v
character(len=maxvarlen), parameter, public var_tb
character(len=maxvarlen), parameter, public var_ts
character(len=maxvarlen), parameter, public var_u
character(len=maxvarlen), parameter, public var_ps
character(len=maxvarlen), parameter, public var_tv