UFO
ufo_gnssro_2d_locs_mod.F90
Go to the documentation of this file.
2 
3 use iso_c_binding
4 use fckit_log_module, only : fckit_log
5 use kinds, only : kind_real
6 use ufo_locs_mod
8 
9 contains
10 
11 !-------------------------------------------------------------------------
12 !-------------------------------------------------------------------------
13 subroutine ufo_gnssro_2d_locs_init(self, locs, obss, t1, t2)
14  use kinds
15  use datetime_mod
16  use obsspace_mod
18 
19  implicit none
20 
21  type(ufo_locs), intent(inout) :: locs
22  class(ufo_gnssro_bndropp2d), intent(inout) :: self
23  type(c_ptr), value, intent(in) :: obss
24  type(datetime), intent(in) :: t1, t2
25 
26  character(len=*),parameter :: myname = "ufo_gnssro_2d_locs_init"
27  integer, parameter :: max_string = 800
28  character(max_string) :: err_msg
29 
30  integer :: i, j, tw_nlocs,nlocs
31  integer, dimension(:), allocatable :: tw_indx
32  real(kind_real), dimension(:), allocatable :: lon, lat
33  type(datetime), dimension(:), allocatable :: date_time
34 
35 ! gnss ro data 2d location
36  real(kind_real), dimension(:), allocatable :: obsazim
37  real(kind_real), dimension(self%roconf%n_horiz) :: plat_2d, plon_2d
38  integer :: kerror, n_horiz
39  real(kind_real) :: dtheta
40 
41  dtheta = self%roconf%dtheta
42  n_horiz = self%roconf%n_horiz
43  ! Local copies pre binning
44  nlocs = obsspace_get_nlocs(obss)
45 
46  allocate(date_time(nlocs), lon(nlocs), lat(nlocs))
47 
48  call obsspace_get_db(obss, "MetaData", "datetime", date_time)
49  call obsspace_get_db(obss, "MetaData", "longitude", lon)
50  call obsspace_get_db(obss, "MetaData", "latitude", lat)
51 
52  ! Generate the timing window indices
53  allocate(tw_indx(nlocs))
54  tw_nlocs = 0
55  do i = 1, nlocs
56  if (date_time(i) > t1 .and. date_time(i) <= t2) then
57  tw_nlocs = tw_nlocs + 1
58  tw_indx(tw_nlocs) = i
59  endif
60  enddo
61 
62  allocate(obsazim(nlocs))
63  if (obsspace_has(obss,"ObsValue","bending_angle")) then
64  if (obsspace_has(obss, "MetaData", "sensor_azimuth_angle")) then
65  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", obsazim)
66  else
67  write(err_msg,*) myname, ' error: sensor_azimuth_angle not found'
68  call abor1_ftn(err_msg)
69  endif
70  endif
71 
72  !Setup ufo 2d locations
73  call ufo_locs_setup(locs, tw_nlocs*n_horiz)
74  do i = 1, tw_nlocs
75  call ropp_fm_2d_plane(lat(tw_indx(i)),lon(tw_indx(i)),obsazim(tw_indx(i)),dtheta,n_horiz,plat_2d,plon_2d,kerror)
76  locs%lon( (i-1)*n_horiz+1 : i*n_horiz) = plon_2d
77  locs%lat( (i-1)*n_horiz+1 : i*n_horiz) = plat_2d
78  ! save ufo_locs to self
79  self%obsLat2d( (i-1)*n_horiz+1 : i*n_horiz) = locs%lat( (i-1)*n_horiz+1 : i*n_horiz)
80  self%obsLon2d( (i-1)*n_horiz+1 : i*n_horiz) = locs%lon( (i-1)*n_horiz+1 : i*n_horiz)
81  do j = 1, n_horiz
82  locs%indx((i-1)*n_horiz+j) = (tw_indx(i)-1)*n_horiz+j
83  locs%time((i-1)*n_horiz+j) = date_time(tw_indx(i))
84  end do
85  end do
86 
87  do i = 1, nlocs
88  call datetime_delete(date_time(i))
89  enddo
90  deallocate(date_time, lon, lat, tw_indx, obsazim)
91 
92 end subroutine ufo_gnssro_2d_locs_init
93 
94 !----------------------------
95 end module ufo_gnssro_2d_locs_mod
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
ufo_gnssro_bndropp2d_mod::ufo_gnssro_bndropp2d
Fortran derived type for gnssro trajectory.
Definition: ufo_gnssro_bndropp2d_mod.F90:33
ufo_gnssro_2d_locs_mod::ufo_gnssro_2d_locs_init
subroutine, public ufo_gnssro_2d_locs_init(self, locs, obss, t1, t2)
Definition: ufo_gnssro_2d_locs_mod.F90:14
ufo_locs_mod::ufo_locs_setup
subroutine, public ufo_locs_setup(self, nlocs)
Definition: ufo_locs_mod.F90:81
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_gnssro_2d_locs_mod
Definition: ufo_gnssro_2d_locs_mod.F90:1
ufo_locs_mod
Fortran module handling observation locations.
Definition: ufo_locs_mod.F90:9
ufo_locs_mod::ufo_locs
Fortran derived type to hold observation locations.
Definition: ufo_locs_mod.F90:25