Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / convertor / decode_l2_airs / module_read_airs.f90
blob3c631e91608ae9fbf523c71489b4c07ce082dd85
2 ! Based on the routine airs_ret_rdr supplied through the AIRS web site
4 module read_airs
6 integer, parameter :: AIRS_RET_GEOXTRACK = 30
7 integer, parameter :: AIRS_RET_GEOTRACK = 45
8 integer, parameter :: AIRS_RET_STDPRESSURELEV = 28
9 integer, parameter :: AIRS_RET_STDPRESSURELAY = 28
10 integer, parameter :: AIRS_RET_H2OPRESSURELEV = 15 ! new in V5
11 integer, parameter :: AIRS_RET_H2OPRESSURELAY = 14 ! new in V5
12 integer, parameter :: AIRS_RET_AIRSXTRACK = 3
13 integer, parameter :: AIRS_RET_AIRSTRACK = 3
14 integer, parameter :: AIRS_RET_CLOUD = 2
15 integer, parameter :: AIRS_RET_CHANAMSUA = 15
16 integer, parameter :: AIRS_RET_CHANHSB = 5
17 integer, parameter :: AIRS_RET_MWHINGESURF = 7
18 integer, parameter :: AIRS_RET_HINGESURF = 100
20 type airs_ret_gran_t
21 double precision :: start_Latitude
22 double precision :: start_Longitude
23 double precision :: start_Time
24 double precision :: end_Latitude
25 double precision :: end_Longitude
26 double precision :: end_Time
27 integer :: start_year
28 integer :: start_month
29 integer :: start_day
30 integer :: start_hour
31 integer :: start_minute
32 real :: start_sec
34 double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Latitude
35 double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Longitude
36 double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Time
38 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nStd_mid_top_bndry
39 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nStd_bot_mid_bndry
40 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: RetQAFlag
41 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Top
42 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Mid
43 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Bot
44 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Cloud_OLR
45 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_H2O
46 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Surf
47 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nBestStd ! new in V5
48 integer*2, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nGoodStd ! new in V5
49 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: PBest ! new in V5
50 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: PGood ! new in V5
52 integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nSurfStd
53 real, dimension(AIRS_RET_STDPRESSURELEV) :: pressStd
54 real, dimension(AIRS_RET_H2OPRESSURELEV) :: pressH2O ! new in V5
55 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: TSurfStd
56 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: TSurfAir
57 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: PSurfStd
58 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Press_bot_mid_bndry
59 real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Press_mid_top_bndry
60 real, dimension(AIRS_RET_STDPRESSURELAY,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: H2OMMRStd
61 real, dimension(AIRS_RET_STDPRESSURELAY,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: H2OMMRStdErr
62 real, dimension(AIRS_RET_STDPRESSURELEV,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: TAirStd
63 real, dimension(AIRS_RET_STDPRESSURELEV,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: GP_Height
65 character (len=256) :: processing_level
66 character (len=256) :: instrument
67 character (len=256) :: DayNightFlag
68 character (len=256) :: AutomaticQAFlag
69 character (len=256) :: node_type
70 character (len=256) :: granules_present
72 integer :: NumTotalData
73 integer :: NumProcessData
74 integer :: NumSpecialData
75 integer :: NumBadData
76 integer :: NumMissingData
77 integer :: NumLandSurface
78 integer :: NumOceanSurface
80 end type airs_ret_gran_t
83 contains
86 subroutine airs_ret_rdr(file_name, airs_ret_gran, version)
88 implicit none
90 ! Arguments
91 character (len=*), intent(in) :: file_name
92 type (airs_ret_gran_t), intent(inout) :: airs_ret_gran
93 character (len=2), intent(out) :: version ! V4 or V5
95 ! Local variables
96 integer :: statn ! HDF-EOS status. 0 for success
97 integer :: fid ! HDF-EOS file ID
98 integer :: swid ! HDF-EOS swath ID
99 integer :: nchar ! Number of characters
100 integer :: nswath ! Number of swaths
101 character (len=256) :: swathname ! Name of swath
102 character (len=256) :: dimnames ! Name of dimensions
103 integer, dimension(10) :: start ! start of each dimensions for Swath I/O
104 ! 0 => start with first element
105 integer, dimension(10) :: stride ! stride of each dimensions for Swath I/O
106 ! 1 => use every element
107 integer, dimension(10) :: edge ! size of each dimension for swath I/O
108 ! will be set for each individual read
109 integer :: swopen, swinqswath, swattach, swinqdims
110 integer :: swrdfld, swrdattr
111 integer :: swdetach, swclose
113 start = 0
114 stride = 1
117 ! Open
119 fid = swopen(file_name, 1)
120 if (fid == -1) then
121 write(6,*) 'Error ', fid, ' opening file ', file_name
122 stop
123 end if
126 ! Get name of swath(s)
128 nswath = swinqswath(file_name, swathname, nchar)
129 if (nswath /= 1) then
130 write(6,*) 'swinqswath found ', nswath, ' swaths for file ', file_name, ' Need exactly 1'
131 stop
132 end if
135 ! There's exactly one swath. Make sure it is the right one.
137 if (swathname /= 'L2_Standard_atmospheric&surface_product') then
138 write(6,*) 'Error: bad swath name ', swathname, ' in file ', file_name
139 write(6,*) 'Expected L2_Standard_atmospheric&surface_product'
140 stop
141 end if
144 ! Attach to (open) the one swath.
146 swid = swattach(fid, swathname)
147 if (swid == -1) then
148 write(6,*) 'Failed to attach to swath ', swathname,' in file ', file_name
149 stop
150 end if
153 ! Read dimension names
155 statn = swinqdims(swid, dimnames, nchar)
156 if ( index(dimnames,'H2OPressureLev') > 0 ) then
157 version = 'V5'
158 else
159 version = 'V4'
160 end if
161 write(6,*) 'Processing version ', version, ' file'
164 ! Read attributes
166 statn = swrdattr(swid, 'processing_level', airs_ret_gran%processing_level)
167 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute processing_level'
169 statn = swrdattr(swid, 'instrument', airs_ret_gran%instrument)
170 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute instrument'
172 statn = swrdattr(swid, 'DayNightFlag', airs_ret_gran%DayNightFlag)
173 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute DayNightFlag'
174 airs_ret_gran%DayNightFlag = airs_ret_gran%DayNightFlag(1:index(airs_ret_gran%DayNightFlag,char(0))-1)
176 statn = swrdattr(swid, 'AutomaticQAFlag', airs_ret_gran%AutomaticQAFlag)
177 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute AutomaticQAFlag'
179 statn = swrdattr(swid, 'NumTotalData', airs_ret_gran%NumTotalData)
180 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumTotalData'
182 statn = swrdattr(swid, 'NumProcessData', airs_ret_gran%NumProcessData)
183 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumProcessData'
185 statn = swrdattr(swid, 'NumSpecialData', airs_ret_gran%NumSpecialData)
186 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumSpecialData'
188 statn = swrdattr(swid, 'NumBadData', airs_ret_gran%NumBadData)
189 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumBadData'
191 statn = swrdattr(swid, 'NumMissingData', airs_ret_gran%NumMissingData)
192 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumMissingData'
194 statn = swrdattr(swid, 'NumLandSurface', airs_ret_gran%NumLandSurface)
195 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumLandSurface'
197 statn = swrdattr(swid, 'NumOceanSurface', airs_ret_gran%NumOceanSurface)
198 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumOceanSurface'
200 statn = swrdattr(swid, 'node_type', airs_ret_gran%node_type)
201 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute node_type'
203 statn = swrdattr(swid, 'start_year', airs_ret_gran%start_year)
204 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_year'
206 statn = swrdattr(swid, 'start_month', airs_ret_gran%start_month)
207 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_month'
209 statn = swrdattr(swid, 'start_day', airs_ret_gran%start_day)
210 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_day'
212 statn = swrdattr(swid, 'start_hour', airs_ret_gran%start_hour)
213 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_hour'
215 statn = swrdattr(swid, 'start_minute', airs_ret_gran%start_minute)
216 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_minute'
218 statn = swrdattr(swid, 'start_sec', airs_ret_gran%start_sec)
219 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_sec'
221 statn = swrdattr(swid, 'start_Latitude', airs_ret_gran%start_Latitude)
222 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Latitude'
224 statn = swrdattr(swid, 'start_Longitude', airs_ret_gran%start_Longitude)
225 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Longitude'
227 statn = swrdattr(swid, 'start_Time', airs_ret_gran%start_Time)
228 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Time'
230 statn = swrdattr(swid, 'end_Latitude', airs_ret_gran%end_Latitude)
231 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Latitude'
233 statn = swrdattr(swid, 'end_Longitude', airs_ret_gran%end_Longitude)
234 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Longitude'
236 statn = swrdattr(swid, 'end_Time', airs_ret_gran%end_Time)
237 if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Time'
239 start = 0
240 stride = 1
243 ! Read geolocation fields
245 edge(1) = AIRS_RET_GEOXTRACK
246 edge(2) = AIRS_RET_GEOTRACK
247 statn = swrdfld(swid, 'Latitude', start, stride, edge, airs_ret_gran%Latitude)
248 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Latitude'
250 statn = swrdfld(swid, 'Longitude', start, stride, edge, airs_ret_gran%Longitude)
251 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Longitude'
253 statn = swrdfld(swid, 'Time', start, stride, edge, airs_ret_gran%Time)
254 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Time'
256 start = 0
257 stride = 1
260 ! Read data Fields
262 edge(2) = 45
263 edge(1) = 30
264 statn = swrdfld(swid, 'RetQAFlag', start, stride, edge, airs_ret_gran%RetQAFlag)
265 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field RetQAFlag'
267 edge(2) = 45
268 edge(1) = 30
269 statn = swrdfld(swid, 'Press_mid_top_bndry', start, stride, edge, airs_ret_gran%Press_mid_top_bndry)
270 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Press_mid_top_bndry'
272 edge(2) = 45
273 edge(1) = 30
274 statn = swrdfld(swid, 'Press_bot_mid_bndry', start, stride, edge, airs_ret_gran%Press_bot_mid_bndry)
275 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Press_bot_mid_bndry'
277 if ( version == 'V5' ) then
278 edge(2) = 45
279 edge(1) = 30
280 statn = swrdfld(swid, 'PBest', start, stride, edge, airs_ret_gran%PBest)
281 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field PBest'
283 edge(2) = 45
284 edge(1) = 30
285 statn = swrdfld(swid, 'PGood', start, stride, edge, airs_ret_gran%PGood)
286 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field PGood'
288 edge(2) = 45
289 edge(1) = 30
290 statn = swrdfld(swid, 'nBestStd', start, stride, edge, airs_ret_gran%nBestStd)
291 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nBestStd'
293 edge(2) = 45
294 edge(1) = 30
295 statn = swrdfld(swid, 'nGoodStd', start, stride, edge, airs_ret_gran%nGoodStd)
296 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nGoodStd'
298 edge(1) = AIRS_RET_H2OPRESSURELEV
299 statn = swrdfld(swid, 'pressH2O', start, stride, edge, airs_ret_gran%pressH2O)
300 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field pressH2O'
302 end if
304 edge(2) = 45
305 edge(1) = 30
306 statn = swrdfld(swid, 'Qual_Cloud_OLR', start, stride, edge, airs_ret_gran%Qual_Cloud_OLR)
307 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Cloud_OLR'
309 edge(2) = 45
310 edge(1) = 30
311 statn = swrdfld(swid, 'Qual_H2O', start, stride, edge, airs_ret_gran%Qual_H2O)
312 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_H2O'
314 edge(2) = 45
315 edge(1) = 30
316 statn = swrdfld(swid, 'Qual_Temp_Profile_Top', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Top)
317 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Top'
319 edge(2) = 45
320 edge(1) = 30
321 statn = swrdfld(swid, 'Qual_Temp_Profile_Mid', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Mid)
322 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Mid'
324 edge(2) = 45
325 edge(1) = 30
326 statn = swrdfld(swid, 'Qual_Temp_Profile_Bot', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Bot)
327 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Bot'
329 edge(2) = 45
330 edge(1) = 30
331 statn = swrdfld(swid, 'Qual_Surf', start, stride, edge, airs_ret_gran%Qual_Surf)
332 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Surf'
334 edge(2) = 45
335 edge(1) = 30
336 statn = swrdfld(swid, 'PSurfStd', start, stride, edge, airs_ret_gran%PSurfStd)
337 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field PSurfStd'
339 edge(2) = 45
340 edge(1) = 30
341 statn = swrdfld(swid, 'nSurfStd', start, stride, edge, airs_ret_gran%nSurfStd)
342 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nSurfStd'
344 edge(1) = 28
345 statn = swrdfld(swid, 'pressStd', start, stride, edge, airs_ret_gran%pressStd)
346 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field pressStd'
348 edge(2) = 45
349 edge(1) = 30
350 statn = swrdfld(swid, 'nStd_mid_top_bndry', start, stride, edge, airs_ret_gran%nStd_mid_top_bndry)
351 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nStd_mid_top_bndry'
353 edge(2) = 45
354 edge(1) = 30
355 statn = swrdfld(swid, 'nStd_bot_mid_bndry', start, stride, edge, airs_ret_gran%nStd_bot_mid_bndry)
356 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nStd_bot_mid_bndry'
358 edge(2) = 45
359 edge(1) = 30
360 statn = swrdfld(swid, 'TSurfStd', start, stride, edge, airs_ret_gran%TSurfStd)
361 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TSurfStd'
363 edge(2) = 45
364 edge(1) = 30
365 statn = swrdfld(swid, 'TSurfAir', start, stride, edge, airs_ret_gran%TSurfAir)
366 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TSurfAir'
368 edge(3) = 45
369 edge(2) = 30
370 edge(1) = 28
371 statn = swrdfld(swid, 'TAirStd', start, stride, edge, airs_ret_gran%TAirStd)
372 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TAirStd'
374 edge(3) = 45
375 edge(2) = 30
376 edge(1) = 28
377 statn = swrdfld(swid, 'GP_Height', start, stride, edge, airs_ret_gran%GP_Height)
378 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field GP_Height'
380 edge(3) = 45
381 edge(2) = 30
382 if ( version == 'V5' ) then
383 edge(1) = AIRS_RET_H2OPRESSURELAY
384 else
385 edge(1) = AIRS_RET_STDPRESSURELAY
386 end if
387 airs_ret_gran%H2OMMRStd = -9999.0 ! initialize
388 statn = swrdfld(swid, 'H2OMMRStd', start, stride, edge, airs_ret_gran%H2OMMRStd)
389 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field H2OMMRStd'
391 edge(3) = 45
392 edge(2) = 30
393 if ( version == 'V5' ) then
394 edge(1) = AIRS_RET_H2OPRESSURELAY
395 else
396 edge(1) = AIRS_RET_STDPRESSURELAY
397 end if
398 airs_ret_gran%H2OMMRStdErr = -9999.0 ! initialize
399 statn = swrdfld(swid, 'H2OMMRStdErr', start, stride, edge, airs_ret_gran%H2OMMRStdErr)
400 if (statn /= 0) write(6,*) 'Error ', statn, ' reading field H2OMMRStdErr'
402 ! Final clean-up
404 statn = swdetach(swid)
405 ! if (statn /= 0) write(6,*) 'Error detaching from input file ', file_name
407 statn = swclose(fid)
408 ! if (statn /= 0) write(6,*) 'Error closing input file ', file_name
410 end subroutine airs_ret_rdr
412 end module read_airs