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