UFO
gnssro_mod_obserror.F90
Go to the documentation of this file.
1 !==========================================================================
3 !==========================================================================
4 
5 use kinds
7 
8 contains
9 subroutine bending_angle_obserr_ecmwf(obsImpH, obsValue, nobs, obsErr, QCflags, missing)
10 implicit none
11 integer, intent(in) :: nobs
12 real(kind_real), dimension(nobs),intent(in) :: obsImpH, obsValue
13 integer(c_int), dimension(nobs),intent(in) :: QCflags(:)
14 real(kind_real), dimension(nobs),intent(out) :: obsErr
15 real(kind_real) :: H_km, missing
16 integer :: i
17 
18 obserr = missing
19 
20 do i = 1, nobs
21 
22 if (qcflags(i) .eq. 0) then
23 
24  h_km = obsimph(i)/1000.0_kind_real
25 
26  if ( h_km <= 10.0 ) then
27  obserr(i) = (h_km*1.25 + (10-h_km)*20)/10.0
28  obserr(i) = obserr(i)/100.0*obsvalue(i)
29  else if ( h_km > 10.0 .and. h_km <= 32.0 ) then
30  obserr(i) = 1.25/100.0*obsvalue(i)
31  else
32  obserr(i) = 3.0*1e-6
33  end if
34 end if
35 
36 end do
37 
38 end subroutine bending_angle_obserr_ecmwf
39 !----------------------------------------
40 
41 subroutine bending_angle_obserr_nrl(obsLat, obsImpH, obsValue, nobs, obsErr, QCflags, missing)
42 implicit none
43 integer, intent(in) :: nobs
44 real(kind_real), dimension(nobs),intent(in) :: obsImpH, obsValue, obsLat
45 integer(c_int), dimension(nobs),intent(in) :: QCflags(:)
46 real(kind_real), dimension(nobs),intent(out) :: obsErr
47 real(kind_real):: H_m, missing
48 real(kind_real):: lat_in_rad,trop_proxy,damping_factor,errfac
49 real(kind_real):: max_sfc_error, min_ba_error
50 
51 integer :: i
52 
53 obserr = missing
54 max_sfc_error = 20.0 ! %
55 min_ba_error = 1.25 ! %
56 
57 do i = 1, nobs
58 
59 if (qcflags(i) .eq. 0) then
60 
61  h_m = obsimph(i)
62  lat_in_rad = deg2rad * obslat(i)
63  trop_proxy = 8666.66 + 3333.33*cos(2.0*lat_in_rad)
64  damping_factor = 0.66 + cos(lat_in_rad)/3.0
65  errfac = max_sfc_error*damping_factor*(trop_proxy - h_m)/trop_proxy
66  errfac = max(min_ba_error, errfac)
67  obserr(i) = max(obsvalue(i)*errfac/100.0, 3.0*1e-6) ! noise floor at top of RO profile
68 
69 end if
70 
71 end do
72 
73 end subroutine bending_angle_obserr_nrl
74 !--------------------------------------
75 
76 subroutine bending_angle_obserr_nbam(obsLat, obsImpH, obsSaid, nobs, obsErr, QCflags, missing)
77 implicit none
78 integer, intent(in) :: nobs
79 real(kind_real), dimension(nobs),intent(in) :: obsImpH, obsLat
80 integer(c_int), dimension(nobs),intent(in) :: obsSaid, QCflags(:)
81 real(kind_real), dimension(nobs),intent(out) :: obsErr
82 real(kind_real) :: H_km, missing
83 
84 integer :: i
85 
86 obserr = missing
87 
88 do i = 1, nobs
89 
90 if (qcflags(i) .eq. 0) then
91 
92  h_km = obsimph(i)/1000.0_kind_real
93  if( (obssaid(i)==41).or.(obssaid(i)==722).or.(obssaid(i)==723).or. &
94  (obssaid(i)==4).or.(obssaid(i)==42).or.(obssaid(i)==3).or. &
95  (obssaid(i)==5).or.(obssaid(i)==821.or.(obssaid(i)==421)).or. &
96  (obssaid(i)==440).or.(obssaid(i)==43)) then
97  if( abs(obslat(i))>= 40.00 ) then
98  if(h_km>12.0) then
99  obserr(i)=0.19032 +0.287535 *h_km-0.00260813*h_km**2
100  else
101  obserr(i)=-3.20978 +1.26964 *h_km-0.0622538 *h_km**2
102  endif
103  else
104  if(h_km>18.) then
105  obserr(i)=-1.87788 +0.354718 *h_km-0.00313189 *h_km**2
106  else
107  obserr(i)=-2.41024 +0.806594 *h_km-0.027257 *h_km**2
108  endif
109  endif
110 
111  else !!!! CDAAC processing
112  if( abs(obslat(i))>= 40.00 ) then
113  if ( h_km>12.00 ) then
114  obserr(i)=-0.685627 +0.377174 *h_km-0.00421934 *h_km**2
115  else
116  obserr(i)=-3.27737 +1.20003 *h_km-0.0558024 *h_km**2
117  endif
118  else
119  if( h_km>18.00 ) then
120  obserr(i)=-2.73867 +0.447663 *h_km-0.00475603 *h_km**2
121  else
122  obserr(i)=-3.45303 +0.908216 *h_km-0.0293331 *h_km**2
123  endif
124  endif
125  obserr(i) = 0.001 /abs(exp(obserr(i)))
126 
127  endif
128 
129 end if
130 
131 end do
132 
133 end subroutine bending_angle_obserr_nbam
134 !---------------------------------------
135 
136 subroutine refractivity_obserr_nbam(obsLat, obsZ, nobs, obsErr, QCflags,missing)
137 implicit none
138 integer, intent(in) :: nobs
139 real(kind_real), dimension(nobs),intent(in) :: obsLat, obsZ
140 real(kind_real), dimension(nobs),intent(out) :: obsErr
141 integer(c_int), dimension(nobs),intent(in) :: QCflags(:)
142 real(kind_real) :: H_km, missing
143 
144 integer :: i
145 
146 obserr = missing
147 
148 do i = 1, nobs
149 
150  if (qcflags(i) .eq. 0) then
151  h_km = obsz(i)/1000.0_kind_real
152  if( abs(obslat(i))>= 20.0 ) then
153  obserr(i)=-1.321+0.341*h_km-0.005*h_km**2
154  else
155  if(h_km > 10.0) then
156  obserr(i)=2.013-0.060*h_km+0.0045*h_km**2
157  else
158  obserr(i)=-1.18+0.058*h_km+0.025*h_km**2
159  endif
160  endif
161  obserr(i) = 1.0_kind_real/abs(exp(obserr(i)))
162  end if
163 end do
164 
165 end subroutine refractivity_obserr_nbam
166 
167 end module gnssro_mod_obserror
168 
gnssro_mod_obserror::refractivity_obserr_nbam
subroutine refractivity_obserr_nbam(obsLat, obsZ, nobs, obsErr, QCflags, missing)
Definition: gnssro_mod_obserror.F90:137
gnssro_mod_obserror
Definition: gnssro_mod_obserror.F90:2
gnssro_mod_obserror::bending_angle_obserr_nbam
subroutine bending_angle_obserr_nbam(obsLat, obsImpH, obsSaid, nobs, obsErr, QCflags, missing)
Definition: gnssro_mod_obserror.F90:77
gnssro_mod_constants
Definition: gnssro_mod_constants.F90:2
gnssro_mod_obserror::bending_angle_obserr_nrl
subroutine bending_angle_obserr_nrl(obsLat, obsImpH, obsValue, nobs, obsErr, QCflags, missing)
Definition: gnssro_mod_obserror.F90:42
gnssro_mod_obserror::bending_angle_obserr_ecmwf
subroutine bending_angle_obserr_ecmwf(obsImpH, obsValue, nobs, obsErr, QCflags, missing)
Definition: gnssro_mod_obserror.F90:10