UFO
ufo_gnssro_bndropp2d_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2018 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
6 !> Fortran module for gnssro bending angle ropp2d forward operator
7 !> following the ROPP (2018 Aug) implementation
8 
10 
11 use fckit_configuration_module, only: fckit_configuration
12 !use iso_c_binding
13 use kinds
14 use ufo_vars_mod
17 use ufo_basis_mod, only: ufo_basis
20 use obsspace_mod
21 use missing_values_mod
24 
26 use ufo_locs_mod
27 use fckit_log_module, only : fckit_log
28 
29 implicit none
30 public :: ufo_gnssro_bndropp2d
31 private
32  !> Fortran derived type for gnssro trajectory
34  type(gnssro_conf) :: roconf
35  real(kind_real), allocatable :: obslon2d(:), obslat2d(:) !2d location
36  contains
37 
38  procedure :: setup => ufo_gnssro_bndropp2d_setup
39  procedure :: simobs => ufo_gnssro_bndropp2d_simobs
40 end type ufo_gnssro_bndropp2d
41 
42 contains
43 
44 ! ------------------------------------------------------------------------------
45 subroutine ufo_gnssro_bndropp2d_setup(self, f_conf, c_size)
46  implicit none
47  class(ufo_gnssro_bndropp2d), intent(inout) :: self
48  type(fckit_configuration), intent(in) :: f_conf
49  integer, intent(in) :: c_size ! 1d obsspace vector length
50 
51  call gnssro_conf_setup(self%roconf,f_conf)
52 
53  allocate(self%obsLon2d(c_size*self%roconf%n_horiz))
54  allocate(self%obsLat2d(c_size*self%roconf%n_horiz))
55 
56  self%obsLon2d = 0.0
57  self%obsLat2d = 0.0
58 
59 end subroutine ufo_gnssro_bndropp2d_setup
60 
61 ! ------------------------------------------------------------------------------
62 subroutine ufo_gnssro_bndropp2d_simobs(self, geovals, hofx, obss)
63  use ropp_fm_types, only: state2dfm, state1dfm
64  use ropp_fm_types, only: obs1dbangle
65  use typesizes, only: wp => eightbytereal
66  use datetimetypes, only: dp
67 
68  implicit none
69  class(ufo_gnssro_bndropp2d), intent(in) :: self
70  type(ufo_geovals), intent(in) :: geovals
71  real(kind_real), intent(inout) :: hofx(:)
72  type(c_ptr), value, intent(in) :: obss
73  real(c_double) :: missing
74 
75  type(state2dfm) :: x
76  type(state1dfm) :: x1d
77  type(obs1dbangle) :: y
78 
79  character(len=*), parameter :: myname_="ufo_gnssro_bndropp2d_simobs"
80  integer, parameter :: max_string = 800
81  character(max_string) :: err_msg
82  integer :: nlev, nlocs, iobs, nvprof
83  integer :: iflip
84  type(ufo_geoval), pointer :: t, q, prs, gph, gph_sfc
85  real(kind_real), allocatable :: obsImpP(:),obsLocR(:),obsGeoid(:),obsAzim(:) !nlocs
86  real(kind_real), allocatable :: obsLat(:),obsLon(:) !nlocs
87  real(kind_real), allocatable :: obsLonnh(:),obsLatnh(:) !n_horiz
88  integer :: n_horiz
89  real(kind_real) :: dtheta
90  real(kind_real) :: ob_time
91  integer, allocatable, dimension(:) :: ichk
92 
93  n_horiz = self%roconf%n_horiz
94  dtheta = self%roconf%dtheta
95 
96  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs: begin"
97  call fckit_log%info(err_msg)
98 ! check if nlocs is consistent in geovals & hofx
99  if (geovals%nlocs /= size(hofx)*n_horiz) then
100  write(err_msg,*) myname_, ' error: 2d nlocs inconsistent! geovals%nlocs, size(hofx), &
101  and n_horiz are', geovals%nlocs, size(hofx), n_horiz
102  call abor1_ftn(err_msg)
103  endif
104 
105 ! get variables from geovals
106  call ufo_geovals_get_var(geovals, var_ts, t) ! temperature
107  call ufo_geovals_get_var(geovals, var_q, q) ! specific humidity
108  call ufo_geovals_get_var(geovals, var_prs, prs) ! pressure
109  call ufo_geovals_get_var(geovals, var_z, gph) ! geopotential height
110  call ufo_geovals_get_var(geovals, var_sfc_geomz, gph_sfc) ! surface geopotential height
111 
112  missing = missing_value(missing)
113  nlev = t%nval ! number of model levels
114  nlocs = obsspace_get_nlocs(obss)
115 
116  iflip = 0
117  if (prs%vals(1,1) .lt. prs%vals(prs%nval,1) ) then
118  iflip = 1
119  write(err_msg,'(a)') ' ufo_gnssro_bndropp2d_simobs:'//new_line('a')// &
120  ' Model vertical height profile is in descending order,'//new_line('a')// &
121  ' but ROPP requires it to be ascending order, need flip'
122  call fckit_log%info(err_msg)
123  end if
124 
125 ! set obs space struture
126  allocate(obslon(nlocs))
127  allocate(obslat(nlocs))
128  allocate(obsimpp(nlocs))
129  allocate(obslocr(nlocs))
130  allocate(obsgeoid(nlocs))
131  allocate(obsazim(nlocs))
132  allocate(obslatnh(n_horiz))
133  allocate(obslonnh(n_horiz))
134 
135  call obsspace_get_db(obss, "MetaData", "longitude", obslon)
136  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
137  call obsspace_get_db(obss, "MetaData", "impact_parameter", obsimpp)
138  call obsspace_get_db(obss, "MetaData", "earth_radius_of_curvature", obslocr)
139  call obsspace_get_db(obss, "MetaData", "geoid_height_above_reference_ellipsoid", obsgeoid)
140  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", obsazim)
141 
142  nvprof = 1 ! no. of bending angles in profile
143  ob_time = 0.0
144  allocate(ichk(nvprof))
145  ichk(:) = 0
146  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs: begin observation loop, nlocs = ", nlocs
147  call fckit_log%info(err_msg)
148 
149 ! loop through the obs
150  obs_loop: do iobs = 1, nlocs
151 
152  if ( ( obsimpp(iobs)-obslocr(iobs)-obsgeoid(iobs) ) <= self%roconf%top_2d .and. &
153  obsazim(iobs) /= missing ) then
154 
155  obslatnh = self%obsLat2d( (iobs-1)*n_horiz+1:iobs*n_horiz )
156  obslonnh = self%obsLon2d( (iobs-1)*n_horiz+1:iobs*n_horiz )
157  call init_ropp_2d_statevec(obslonnh, obslatnh, &
158  t%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
159  q%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
160  prs%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
161  gph%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
162  nlev,x,n_horiz,dtheta,iflip)
163 
164  call init_ropp_2d_obvec(nvprof, &
165  obsimpp(iobs), &
166  obslat(iobs), &
167  obslon(iobs), &
168  obslocr(iobs), &
169  obsgeoid(iobs), &
170  y)
171 
172  call ropp_fm_bangle_2d(x,y)
173 
174  else ! apply ropp1d above top_2d
175 
176  call init_ropp_1d_statevec(ob_time, &
177  obslon(iobs), &
178  obslat(iobs), &
179  t%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
180  q%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
181  prs%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
182  gph%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
183  nlev, &
184  gph_sfc%vals(1,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
185  x1d, iflip)
186 
187  call init_ropp_1d_obvec(nvprof, &
188  obsimpp(iobs), &
189  ichk, ob_time, &
190  obslat(iobs), &
191  obslon(iobs), &
192  obslocr(iobs), &
193  obsgeoid(iobs), &
194  y)
195 
196  call ropp_fm_bangle_1d(x1d,y)
197 
198  end if
199 
200 ! hack -- handling ropp missing value
201  if (y%bangle(nvprof) .lt. -900.0_wp ) then
202  hofx(iobs) = missing
203  y%bangle(nvprof) = missing
204  else
205  hofx(iobs) = y%bangle(nvprof) ! nvprof is just one point
206  end if
207 ! hack -- handling ropp missing value
208 
209 ! deallocate ropp structures
210  if ( ( obsimpp(iobs)-obslocr(iobs)-obsgeoid(iobs) ) <= self%roconf%top_2d ) then
211  call ropp_tidy_up_2d(x,y)
212  else
213  call ropp_tidy_up_1d(x1d,y)
214  end if
215 
216  end do obs_loop
217 
218  deallocate(obslat)
219  deallocate(obslon)
220  deallocate(obsimpp)
221  deallocate(obslocr)
222  deallocate(obsgeoid)
223  deallocate(obsazim)
224  deallocate(obslatnh)
225  deallocate(obslonnh)
226  deallocate(ichk)
227 
228  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs: completed"
229  call fckit_log%info(err_msg)
230 
231  return
232 
233 end subroutine ufo_gnssro_bndropp2d_simobs
234 ! ------------------------------------------------------------------------------
235 
236 end module ufo_gnssro_bndropp2d_mod
lag_interp_mod
Fortran module to prepare for Lagrange polynomial interpolation. based on GSI: lagmod....
Definition: lag_interp.F90:4
ufo_gnssro_ropp1d_utils_mod::init_ropp_1d_statevec
subroutine, public init_ropp_1d_statevec(step_time, rlon, rlat, temp, shum, pres, phi, lm, phi_sfc, x, iflip)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:37
ufo_basis_mod
Definition: ufo_basis_mod.F90:6
ufo_gnssro_ropp2d_utils_mod::init_ropp_2d_obvec
subroutine, public init_ropp_2d_obvec(nvprof, obs_impact, rlat, rlon, roc, undulat, y)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:204
ufo_gnssro_bndropp2d_mod::ufo_gnssro_bndropp2d
Fortran derived type for gnssro trajectory.
Definition: ufo_gnssro_bndropp2d_mod.F90:33
ufo_gnssro_bndropp2d_mod::ufo_gnssro_bndropp2d_setup
subroutine ufo_gnssro_bndropp2d_setup(self, f_conf, c_size)
Definition: ufo_gnssro_bndropp2d_mod.F90:46
ufo_gnssro_ropp2d_utils_mod::init_ropp_2d_statevec
subroutine, public init_ropp_2d_statevec(rlon, rlat, temp, shum, pres, phi, lm, x, n_horiz, dtheta, iflip)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:36
lag_interp_mod::lag_interp_smthweights
subroutine, public lag_interp_smthweights(x, xt, aq, bq, w, dw, n)
Definition: lag_interp.F90:174
gnssro_mod_conf::gnssro_conf
Definition: gnssro_mod_conf.F90:14
ufo_gnssro_ropp2d_utils_mod
Fortran module to handle gnssro bending angle observations following the ROPP (2018 Aug) implementati...
Definition: ufo_gnssro_ropp2d_utils_mod.F90:9
vert_interp_mod
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_gnssro_ropp1d_utils_mod::ropp_tidy_up_1d
subroutine, public ropp_tidy_up_1d(x, y)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:382
ufo_geovals_mod_c
Definition: GeoVaLs.interface.F90:7
ufo_gnssro_ropp1d_utils_mod::init_ropp_1d_obvec
subroutine, public init_ropp_1d_obvec(nvprof, obs_impact, ichk, ob_time, rlat, rlon, roc, undulat, y)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:237
ufo_vars_mod::var_sfc_geomz
character(len=maxvarlen), parameter, public var_sfc_geomz
Definition: ufo_variables_mod.F90:70
gnssro_mod_conf::gnssro_conf_setup
subroutine, public gnssro_conf_setup(roconf, f_conf)
Definition: gnssro_mod_conf.F90:35
ufo_gnssro_ropp2d_utils_mod::ropp_tidy_up_2d
subroutine, public ropp_tidy_up_2d(x, y)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:351
ufo_basis_mod::ufo_basis
Definition: ufo_basis_mod.F90:12
lag_interp_mod::lag_interp_const
subroutine, public lag_interp_const(q, x, n)
Definition: lag_interp.F90:24
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_gnssro_bndropp2d_mod::ufo_gnssro_bndropp2d_simobs
subroutine ufo_gnssro_bndropp2d_simobs(self, geovals, hofx, obss)
Definition: ufo_gnssro_bndropp2d_mod.F90:63
ufo_vars_mod::var_z
character(len=maxvarlen), parameter, public var_z
Definition: ufo_variables_mod.F90:29
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_gnssro_bndropp2d_mod
Fortran module for gnssro bending angle ropp2d forward operator following the ROPP (2018 Aug) impleme...
Definition: ufo_gnssro_bndropp2d_mod.F90:9
ufo_vars_mod::var_q
character(len=maxvarlen), parameter, public var_q
Definition: ufo_variables_mod.F90:22
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_locs_mod
Fortran module handling observation locations.
Definition: ufo_locs_mod.F90:9
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
gnssro_mod_conf
Definition: gnssro_mod_conf.F90:2
ufo_gnssro_ropp1d_utils_mod
Fortran module to handle gnssro bending angle observations following the ROPP (2018 Aug) implementati...
Definition: ufo_gnssro_ropp1d_utils_mod.F90:9
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25
ufo_geovals_mod_c::ufo_geovals_registry
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
Definition: GeoVaLs.interface.F90:30