module_read_airs.f90

References to this file elsewhere.
1 !
2 ! Based on the routine airs_ret_rdr supplied through the AIRS web site
3 !
4 module read_airs
5 
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_AIRSXTRACK     =   3
11    integer, parameter :: AIRS_RET_AIRSTRACK      =   3
12    integer, parameter :: AIRS_RET_CLOUD          =   2
13    integer, parameter :: AIRS_RET_CHANAMSUA      =  15
14    integer, parameter :: AIRS_RET_CHANHSB        =   5
15    integer, parameter :: AIRS_RET_MWHINGESURF    =   7
16    integer, parameter :: AIRS_RET_HINGESURF      = 100
17 
18    type airs_ret_gran_t
19       double precision :: start_Latitude
20       double precision :: start_Longitude
21       double precision :: start_Time
22       double precision :: end_Latitude
23       double precision :: end_Longitude
24       double precision :: end_Time
25       integer :: start_year
26       integer :: start_month
27       integer :: start_day
28       integer :: start_hour
29       integer :: start_minute
30       real :: start_sec
31 
32       double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Latitude
33       double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Longitude
34       double precision, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Time
35 
36       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nStd_mid_top_bndry
37       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nStd_bot_mid_bndry
38       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: RetQAFlag
39       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Top
40       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Mid
41       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Temp_Profile_Bot
42       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Cloud_OLR
43       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_H2O
44       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: Qual_Surf
45 
46       integer, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: nSurfStd
47       real, dimension(AIRS_RET_STDPRESSURELEV)                 :: pressStd
48       real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK)    :: TSurfStd
49       real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK)    :: TSurfAir
50       real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK)    :: PSurfStd
51       real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK)    :: Press_bot_mid_bndry
52       real, dimension(AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK)    :: Press_mid_top_bndry
53       real, dimension(AIRS_RET_STDPRESSURELAY,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: H2OMMRStd
54       real, dimension(AIRS_RET_STDPRESSURELEV,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: TAirStd
55       real, dimension(AIRS_RET_STDPRESSURELEV,AIRS_RET_GEOXTRACK,AIRS_RET_GEOTRACK) :: GP_Height
56 
57       character (len=256) :: processing_level
58       character (len=256) :: instrument
59       character (len=256) :: DayNightFlag
60       character (len=256) :: AutomaticQAFlag
61       character (len=256) :: node_type
62       character (len=256) :: granules_present
63 
64       integer :: NumTotalData
65       integer :: NumProcessData
66       integer :: NumSpecialData
67       integer :: NumBadData
68       integer :: NumMissingData
69       integer :: NumLandSurface
70       integer :: NumOceanSurface
71 
72    end type airs_ret_gran_t
73 
74 
75    contains
76 
77 
78    subroutine airs_ret_rdr(file_name, airs_ret_gran)
79 
80       implicit none
81 
82       ! Arguments
83       character (len=*) :: file_name
84       type (airs_ret_gran_t) :: airs_ret_gran
85    
86       ! Local variables
87       integer :: statn           ! HDF-EOS status. 0 for success
88       integer :: fid             ! HDF-EOS file ID
89       integer :: swid            ! HDF-EOS swath ID
90       integer :: nchar           ! Number of characters
91       integer :: nswath          ! Number of swaths
92       character (len=256) :: swathname   ! Name of swath
93       integer, dimension(10) :: start  ! start of each dimensions for Swath I/O
94                                        ! 0 => start with first element
95       integer, dimension(10) :: stride ! stride of each dimensions for Swath I/O
96                                        ! 1 => use every element
97       integer, dimension(10) :: edge   ! size of each dimension for swath I/O
98                                        ! will be set for each individual read
99       integer :: swopen, swinqswath, swattach
100       integer :: swrdfld, swrdattr
101       integer :: swdetach, swclose
102 
103       start = 0
104       stride = 1
105    
106       !
107       ! Open
108       !
109       fid = swopen(file_name, 1)
110       if (fid == -1) then
111          write(6,*) 'Error ', fid, ' opening file ', file_name
112          stop
113       end if
114    
115       !
116       ! Get name of swath(s)
117       !
118       nswath = swinqswath(file_name, swathname, nchar)
119       if (nswath /= 1) then
120          write(6,*) 'swinqswath found ', nswath, ' swaths for file ', file_name, ' Need exactly 1'
121          stop
122       end if
123    
124       !
125       ! There's exactly one swath.  Make sure it is the right one.
126       !
127       if (swathname /= 'L2_Standard_atmospheric&surface_product') then
128          write(6,*) 'Error: bad swath name ', swathname, ' in file ', file_name
129          write(6,*) 'Expected L2_Standard_atmospheric&surface_product'
130          stop
131       end if
132    
133       !
134       ! Attach to (open) the one swath.
135       !
136       swid = swattach(fid, swathname)
137       if (swid == -1) then
138         write(6,*) 'Failed to attach to swath ', swathname,' in file ', file_name
139         stop
140       end if
141    
142       !
143       ! Read attributes
144       !
145       statn = swrdattr(swid, 'processing_level', airs_ret_gran%processing_level)
146       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute processing_level'
147    
148       statn = swrdattr(swid, 'instrument', airs_ret_gran%instrument)
149       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute instrument'
150    
151       statn = swrdattr(swid, 'DayNightFlag', airs_ret_gran%DayNightFlag)
152       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute DayNightFlag'
153    
154       statn = swrdattr(swid, 'AutomaticQAFlag', airs_ret_gran%AutomaticQAFlag)
155       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute AutomaticQAFlag'
156    
157       statn = swrdattr(swid, 'NumTotalData', airs_ret_gran%NumTotalData)
158       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumTotalData'
159    
160       statn = swrdattr(swid, 'NumProcessData', airs_ret_gran%NumProcessData)
161       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumProcessData'
162    
163       statn = swrdattr(swid, 'NumSpecialData', airs_ret_gran%NumSpecialData)
164       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumSpecialData'
165    
166       statn = swrdattr(swid, 'NumBadData', airs_ret_gran%NumBadData)
167       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumBadData'
168    
169       statn = swrdattr(swid, 'NumMissingData', airs_ret_gran%NumMissingData)
170       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumMissingData'
171    
172       statn = swrdattr(swid, 'NumLandSurface', airs_ret_gran%NumLandSurface)
173       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumLandSurface'
174    
175       statn = swrdattr(swid, 'NumOceanSurface', airs_ret_gran%NumOceanSurface)
176       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute NumOceanSurface'
177    
178       statn = swrdattr(swid, 'node_type', airs_ret_gran%node_type)
179       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute node_type'
180    
181       statn = swrdattr(swid, 'start_year', airs_ret_gran%start_year)
182       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_year'
183    
184       statn = swrdattr(swid, 'start_month', airs_ret_gran%start_month)
185       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_month'
186    
187       statn = swrdattr(swid, 'start_day', airs_ret_gran%start_day)
188       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_day'
189    
190       statn = swrdattr(swid, 'start_hour', airs_ret_gran%start_hour)
191       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_hour'
192    
193       statn = swrdattr(swid, 'start_minute', airs_ret_gran%start_minute)
194       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_minute'
195    
196       statn = swrdattr(swid, 'start_sec', airs_ret_gran%start_sec)
197       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_sec'
198    
199       statn = swrdattr(swid, 'start_Latitude', airs_ret_gran%start_Latitude)
200       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Latitude'
201    
202       statn = swrdattr(swid, 'start_Longitude', airs_ret_gran%start_Longitude)
203       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Longitude'
204    
205       statn = swrdattr(swid, 'start_Time', airs_ret_gran%start_Time)
206       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute start_Time'
207    
208       statn = swrdattr(swid, 'end_Latitude', airs_ret_gran%end_Latitude)
209       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Latitude'
210    
211       statn = swrdattr(swid, 'end_Longitude', airs_ret_gran%end_Longitude)
212       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Longitude'
213    
214       statn = swrdattr(swid, 'end_Time', airs_ret_gran%end_Time)
215       if (statn /= 0) write(6,*) 'Error ', statn, ' reading attribute end_Time'
216    
217 
218       !
219       ! Read geolocation fields
220       !
221       edge(1) = AIRS_RET_GEOXTRACK
222       edge(2) = AIRS_RET_GEOTRACK
223       statn = swrdfld(swid, 'Latitude', start, stride, edge, airs_ret_gran%Latitude)
224       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Latitude'
225    
226       statn = swrdfld(swid, 'Longitude', start, stride, edge, airs_ret_gran%Longitude)
227       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Longitude'
228    
229       statn = swrdfld(swid, 'Time', start, stride, edge, airs_ret_gran%Time)
230       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Time'
231    
232    
233       !
234       ! Read data Fields
235       !
236       edge(2) = 45
237       edge(1) = 30
238       statn = swrdfld(swid, 'RetQAFlag', start, stride, edge, airs_ret_gran%RetQAFlag)
239       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field RetQAFlag'
240 
241       edge(2) = 45
242       edge(1) = 30
243       statn = swrdfld(swid, 'Press_mid_top_bndry', start, stride, edge, airs_ret_gran%Press_mid_top_bndry)
244       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Press_mid_top_bndry'
245 
246       edge(2) = 45
247       edge(1) = 30
248       statn = swrdfld(swid, 'Press_bot_mid_bndry', start, stride, edge, airs_ret_gran%Press_bot_mid_bndry)
249       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Press_bot_mid_bndry'
250 
251       edge(2) = 45
252       edge(1) = 30
253       statn = swrdfld(swid, 'Qual_Cloud_OLR', start, stride, edge, airs_ret_gran%Qual_Cloud_OLR)
254       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Cloud_OLR'
255    
256       edge(2) = 45
257       edge(1) = 30
258       statn = swrdfld(swid, 'Qual_H2O', start, stride, edge, airs_ret_gran%Qual_H2O)
259       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_H2O'
260    
261       edge(2) = 45
262       edge(1) = 30
263       statn = swrdfld(swid, 'Qual_Temp_Profile_Top', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Top)
264       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Top'
265    
266       edge(2) = 45
267       edge(1) = 30
268       statn = swrdfld(swid, 'Qual_Temp_Profile_Mid', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Mid)
269       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Mid'
270    
271       edge(2) = 45
272       edge(1) = 30
273       statn = swrdfld(swid, 'Qual_Temp_Profile_Bot', start, stride, edge, airs_ret_gran%Qual_Temp_Profile_Bot)
274       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Temp_Profile_Bot'
275    
276       edge(2) = 45
277       edge(1) = 30
278       statn = swrdfld(swid, 'Qual_Surf', start, stride, edge, airs_ret_gran%Qual_Surf)
279       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field Qual_Surf'
280    
281       edge(2) = 45
282       edge(1) = 30
283       statn = swrdfld(swid, 'PSurfStd', start, stride, edge, airs_ret_gran%PSurfStd)
284       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field PSurfStd'
285    
286       edge(2) = 45
287       edge(1) = 30
288       statn = swrdfld(swid, 'nSurfStd', start, stride, edge, airs_ret_gran%nSurfStd)
289       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nSurfStd'
290 
291       edge(1) = 28
292       statn = swrdfld(swid, 'pressStd', start, stride, edge, airs_ret_gran%pressStd)
293       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field pressStd'
294    
295       edge(2) = 45
296       edge(1) = 30
297       statn = swrdfld(swid, 'nStd_mid_top_bndry', start, stride, edge, airs_ret_gran%nStd_mid_top_bndry)
298       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nStd_mid_top_bndry'
299    
300       edge(2) = 45
301       edge(1) = 30
302       statn = swrdfld(swid, 'nStd_bot_mid_bndry', start, stride, edge, airs_ret_gran%nStd_bot_mid_bndry)
303       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field nStd_bot_mid_bndry'
304    
305       edge(2) = 45
306       edge(1) = 30
307       statn = swrdfld(swid, 'TSurfStd', start, stride, edge, airs_ret_gran%TSurfStd)
308       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TSurfStd'
309    
310       edge(2) = 45
311       edge(1) = 30
312       statn = swrdfld(swid, 'TSurfAir', start, stride, edge, airs_ret_gran%TSurfAir)
313       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TSurfAir'
314    
315       edge(3) = 45
316       edge(2) = 30
317       edge(1) = 28
318       statn = swrdfld(swid, 'TAirStd', start, stride, edge, airs_ret_gran%TAirStd)
319       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field TAirStd'
320 
321       edge(3) = 45
322       edge(2) = 30
323       edge(1) = 28
324       statn = swrdfld(swid, 'GP_Height', start, stride, edge, airs_ret_gran%GP_Height)
325       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field GP_Height'
326    
327       edge(3) = 45
328       edge(2) = 30
329       edge(1) = 28
330       statn = swrdfld(swid, 'H2OMMRStd', start, stride, edge, airs_ret_gran%H2OMMRStd)
331       if (statn /= 0) write(6,*) 'Error ', statn, ' reading field H2OMMRStd'
332    
333       !
334       ! Final clean-up
335       !
336       statn = swdetach(swid)
337       if (statn /= 0) write(6,*) 'Error detaching from input file ', file_name
338 
339       statn = swclose(fid)
340       if (statn /= 0) write(6,*) 'Error closing input file ', file_name
341    
342    end subroutine airs_ret_rdr
343    
344 end module read_airs