da_read_ssmi_info.inc

References to this file elsewhere.
1 subroutine da_read_obs_ssmi_info (iunit, ob, xb, xbx)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Read the header of a SSMI GTS observation file
5    !---------------------------------------------------------------------------
6 
7    implicit none
8 
9    integer,        intent (in)  :: iunit
10    type (xb_type), intent (in)  :: xb
11    type (xbx_type),intent (in)  :: xbx     ! Header & non-gridded vars.
12    type (iv_type), intent (out) :: ob
13 
14 
15    integer                      :: iost, i, ii
16    character (LEN = 120)        :: char_ned
17    logical                      :: connected
18 
19    integer                      :: nssmis,nothers
20    integer                      :: ixc, jxc, iproj, idd, maxnes
21    integer                      :: nestix(10) , nestjx(10) , numnc(10) , nesti(10) , nestj(10) 
22    real                         :: phic   , xlonc  , &
23                                    truelat1_3dv, truelat2_3dv, &
24                                    local_ts0    , local_ps0    , local_tlp     , ptop
25    real                         :: dis(10)
26 
27    logical                      :: check_wrong, check_subdomain
28 
29    if (trace_use) call da_trace_entry("da_read_obs_ssmi_info")
30 
31    iost = -99999
32 
33    ! 1. open file
34    !    ---------
35 
36    if (use_ssmiretrievalobs .or. use_ssmitbobs .or. use_ssmt1obs .or. use_ssmt2obs) then
37       open (unit   = iunit,     &
38          FORM   = 'FORMATTED',  &
39          ACCESS = 'SEQUENTIAL', &
40          iostat =  iost,     &
41          STATUS = 'OLD')
42 
43       if (iost /= 0) then
44          call da_warning(__FILE__,__LINE__, (/"Cannot open SSMI file"/))
45          use_ssmiretrievalobs = .false.
46          use_ssmitbobs  = .false.
47          use_ssmt1obs   = .false.
48          use_ssmt2obs   = .false.
49          return
50       end if
51    else
52       return
53    end if
54 
55    rewind (unit = iunit)
56 
57 
58    ! 2.  read header
59    ! ===============
60 
61    ! 2.1 read first line
62    !     ---------------
63 
64    read (unit = iunit, fmt = '(a)', iostat = iost) char_ned
65    if (iost /= 0) then
66       use_ssmiretrievalobs = .false.
67       use_ssmitbobs  = .false.
68       use_ssmt1obs   = .false.
69       use_ssmt2obs   = .false.
70       call da_error(__FILE__,__LINE__, (/"Cannot read SSMI file"/))
71    end if
72 
73    ! 2.3 read NUMBER OF REPORTS
74    !     ---------------------
75 
76    do
77       do i = 0, 120-14
78          select case (char_ned (I+1:I+5))
79          ! Number of observations
80          case ('SSMI ') ; 
81             if (use_ssmiretrievalobs) &
82                read (char_ned (I+9:I+14),'(I6)', iostat = iost) &
83                   ob%nlocal(ssmi_rv)
84             if (use_ssmitbobs) then
85                read (char_ned (I+9:I+14),'(I6)', iostat = iost) ob%nlocal(ssmi_tb)
86             end if
87          case ('OTHER') ; 
88             read (char_ned (I+9:I+14),'(I6)', iostat = iost) nothers
89 
90             ! Geographic area and reference atmosphere definition
91 
92          case ('MISS.') ; 
93              read (char_ned (I+8:I+15),'(F8.0)', iostat = iost) ob % missing
94          case ('PHIC ') ; 
95              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) phic
96          case ('XLONC') ; 
97              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) xlonc
98          case ('true1') ; 
99              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) truelat1_3dv
100          case ('true2') ; 
101              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) truelat2_3dv
102          case ('TS0  ') ; 
103              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) local_ts0
104          case ('TLP  ') ; 
105              read (char_ned (I+8:I+14),'(F7.2)', iostat = iost) local_tlp
106          case ('PTOP ') ; 
107              read (char_ned (I+8:I+14),'(F7.0)', iostat = iost) ptop
108          case ('PS0  ') ; 
109              read (char_ned (I+8:I+14),'(F7.0)', iostat = iost) local_ps0
110          case ('IXC  ') ; 
111              read (char_ned (I+8:I+14),'(I7)', iostat = iost) ixc
112          case ('JXC  ') ; 
113              read (char_ned (I+8:I+14),'(I7)', iostat = iost) jxc
114          case ('IPROJ') ; 
115              read (char_ned (I+8:I+14),'(I7)', iostat = iost) iproj
116          case ('IDD  ') ; 
117              read (char_ned (I+8:I+14),'(I7)', iostat = iost) idd
118          case ('MAXNE') ; 
119              read (char_ned (I+8:I+14),'(I7)', iostat = iost) maxnes
120          case default ; read (char_ned (I+9:I+14),'(I6)', iostat = iost) nssmis
121          end select
122       end do
123 
124       read (unit = iunit, fmt = '(A)', iostat = iost) char_ned
125 
126       if (iost /= 0) then
127          use_ssmiretrievalobs = .false.
128          use_ssmitbobs  = .false.
129          use_ssmt1obs   = .false.
130          use_ssmt2obs   = .false.
131          call da_warning(__FILE__,__LINE__, &
132             (/"Cannot read SSMI file"/))
133          return
134       end if
135       if (char_ned(1:6) == 'NESTIX') exit
136    end do
137 
138    do
139       select case (char_ned (1:6))
140       ! Model domains definition
141 
142       case ('NESTIX') ;
143          call da_read_obs_ssmi_integer_array(nestix, maxnes, 8, 9)
144       case ('NESTJX') ; 
145          call da_read_obs_ssmi_integer_array(nestjx, maxnes, 8, 9)
146       case ('NUMC  ') ; 
147          call da_read_obs_ssmi_integer_array(numnc , maxnes, 8, 9)
148       case ('DIS   ') ; 
149          call da_read_obs_ssmi_real_array   (dis   , maxnes, 8, 9)
150       case ('NESTI ') ; 
151          call da_read_obs_ssmi_integer_array(nesti , maxnes, 8, 9)
152       case ('NESTJ ') ; 
153          call da_read_obs_ssmi_integer_array(nestj , maxnes, 8, 9)
154       end select
155 
156       read (unit = iunit, fmt = '(A)', iostat = iost) char_ned
157       if (char_ned(1:6) == 'INFO  ') exit
158    end do
159 
160    read (unit = iunit, fmt = '(A)', iostat = iost) char_ned
161 
162    if (trace_use) call da_trace_exit("da_read_obs_ssmi_info")
163 
164 contains
165 
166 #include "da_read_obs_ssmi_integer_array.inc"
167 #include "da_read_obs_ssmi_real_array.inc"
168 
169 end subroutine da_read_obs_ssmi_info
170 
171